LDMX Software
EventReadoutProducer.cxx
2
3#include <iostream>
4
5#include "Framework/Exception/Exception.h"
6
7namespace trigscint {
8
9EventReadoutProducer::EventReadoutProducer(const std::string &name,
10 framework::Process &process)
11 : Producer(name, process) {}
12
13void EventReadoutProducer::configure(
15 // Configure this instance of the producer
16 input_collection_ = parameters.get<std::string>("input_collection");
17 input_pass_name_ = parameters.get<std::string>("input_pass_name");
18 output_collection_ = parameters.get<std::string>("output_collection");
19 n_ped_samples_ = parameters.get<int>("number_pedestal_samples");
20 time_shift_ = parameters.get<int>("time_shift");
21 fiber_to_shift_ = parameters.get<int>("fiber_to_shift");
22 verbose_ = parameters.get<bool>("verbose");
23
24 ldmx_log(debug) << "In configure, got parameters:" << "\noutput_collection = "
25 << output_collection_
26 << "\ninput_collection = " << input_collection_
27 << "\ninput_pass_name = " << input_pass_name_
28 << "\nnumber_pedestal_samples = " << n_ped_samples_
29 << "\ntime_shift = " << time_shift_
30 << "\nfiber_to_shift = " << fiber_to_shift_
31 << "\nverbose = " << verbose_;
32}
33
34void EventReadoutProducer::produce(framework::Event &event) {
35 // initialize QIE object for linearizing ADCs
36 SimQIE qie;
37
38 const auto digis{event.getCollection<trigscint::TrigScintQIEDigis>(
39 input_collection_, input_pass_name_)};
40
41 std::vector<trigscint::EventReadout> channel_readout_events;
42 for (const auto &digi : digis) {
44 auto adc{digi.getADC()};
45 auto tdc{digi.getTDC()};
46
47 // copy over from qie digi for convenience
48 out_event.setChanID(digi.getChanID());
49 out_event.setElecID(digi.getElecID());
50 out_event.setTimeSinceSpill(digi.getTimeSinceSpill());
51 // elecID increases monotonically with 8 channels per fiber
52 out_event.setFiberNb(digi.getElecID() / 8);
53 if (out_event.getFiberNb() == fiber_to_shift_)
54 out_event.setTimeOffset(time_shift_);
55
56 out_event.setADC(adc);
57 out_event.setTDC(tdc);
58 std::vector<float> charge;
59 std::vector<float> charge_err;
60
61 float avg_q = 0;
62 float tot_pos_q = 0;
63 int i_s = 0;
64 [[maybe_unused]] int n_pos = 0;
65 float early_ped = 0;
66 for (auto &val : adc) {
67 float q = qie.adc2Q(val);
68 charge.push_back(q);
69 charge_err.push_back(qie.qErr(q));
70 avg_q += q; // charge.back();
71 if (q > 0) {
72 tot_pos_q += q;
73 n_pos++;
74 }
75 if (verbose_)
76 ldmx_log(debug) << "got adc value " << val << " and charge "
77 << q; // qie.ADC2Q(val);
78 if (i_s < n_ped_samples_) early_ped += q;
79 i_s++;
80 }
81 out_event.setQ(charge); // set in proper order before sorting
82 out_event.setQError(charge_err); // set in proper order before sorting
83 early_ped /= n_ped_samples_;
84 out_event.setEarlyPedestal(early_ped);
85
86 // oscillation check. for this the pulse charge needs to be in order, so set
87 // this up now the period is 4. ansatz: an oscillation is a repeated shape.
88 // normalize to maxq=1, in every interval of 4 samples if after
89 // normalization the same numbers are repeated, it's an oscillation edge
90 // case: multiple single PE peaks with that repetition. we probably don't
91 // need to keep those anyway
92 std::vector<float> charge_check = {NULL};
93 float min_charge = 10;
94 // no point in looking at oscillations just around the pedestal.
95 // nearest edge is 10.35 fC //ADC=0 corresponds to -16 fC.
96 float ped = 0;
97 int ped_length = (int)charge.size() / 5;
98 // actually pulse can be up to 12 samples = 2/5*30
99 int ped_offset = 2;
100 // but still skip the lowest (and highest) few
101 // for (int i = pedLength; i < 3*pedLength ; i++) {
103 if (charge.size() > 8) {
104 if (verbose_) ldmx_log(debug) << "going into oscillations check ";
105 for (int i = 3; i < charge.size() - 4; i++) {
106 float max_samp = min_charge;
107 for (int i_q = 0; i_q < 4;
108 i_q++) { // find the local max in the 4 samples
109 // ldmx_log(debug) << "got charge " << charge[i+iQ];
110 if (charge[i + i_q] > max_samp) max_samp = charge[i + i_q];
111 }
112 if (verbose_) ldmx_log(debug) << "got max charge " << max_samp;
113 for (int i_q = 0; i_q < 4; i_q++) // store the locally normalized
114 // numbers. even if the period is
115 // 5 this should work for a while
116 charge_check.push_back(charge[i + i_q] / max_samp);
117 i += 3; // to increment by 4, do 3 here and 1 in the loop
118 }
119 }
120 if (ped_length > 4) {
121 // now calculate the pedestal as the average of the middle half of the
122 // sorted vector
123 std::sort(charge.begin(), charge.end());
124 for (int i = ped_offset; i < 2 * ped_length + ped_offset;
125 i++) { // use 1st and 2nd 5th
126 ped += charge[i];
127 }
128 ped /= 2 * ped_length;
129 }
130 // median: technically only true for odd number of elements but good enough
131 float med_q = charge[(int)charge.size() / 2];
132 float min_q = charge[0];
133 float max_q = charge[charge.size() - 1];
134
135 out_event.setTotQ(tot_pos_q);
136 //-nPos*ped); //store (event) ped subtracted
137 // total charge, before dividing by N -->
138 // actually, ped subtraction makes it confusing
139 // outEvent.setTotQ(totPosQ-adc.size()*ped); //store (event) ped subtracted
140 // total charge, before dividing by N
141 // outEvent.setTotQ(avgQ-adc.size()*ped);
143 avg_q /= adc.size();
144 out_event.setPedestal(ped);
145 out_event.setAvgQ(avg_q);
146 out_event.setMedQ(med_q);
147 out_event.setMinQ(min_q);
148 out_event.setMaxQ(max_q);
149
150 // and the noise as the RMSE of that, same interval as pedestal
151 float diff_sq = 0;
152 // for (int i = pedLength; i < 3*pedLength ; i++) {
153 if (charge.size() > 8) {
154 for (int i = ped_offset; i < 2 * ped_length + ped_offset; i++) {
155 diff_sq += (charge[i] - ped) * (charge[i] - ped);
156 }
157 diff_sq /= 2 * ped_length; // adc.size();
158 }
159 out_event.setNoise(sqrt(diff_sq));
160
161 // oscillation check
162 uint flag_oscillation = 0;
163 if (charge.size() > 8) {
164 // no need to run tedious oscillation check for all-neg channels
165 if (max_q > min_charge) {
166 int max_id = 0;
167 // find the first occurence of a local max
168 for (int i = 0; i < charge_check.size() - 4; i++) {
169 if (charge_check[i] ==
170 1) { // an actual local max has been normalised by its own value
171 max_id = i;
172 // ldmx_log(debug) << "storing max index " <<maxID <<"
173 // and size of vector is " << chargeCheck.size()-4;
174 break;
175 }
176 }
177 int last_match_sample = 0;
178 // start from local max
179 bool do_break = false;
180 for (int i = max_id; i < charge_check.size() - 4; i++) {
181 if (verbose_)
182 ldmx_log(debug)
183 << "Checking how many matching groups of four we can "
184 "find, starting at index "
185 << i;
186
187 for (int i_q = 0; i_q < 4; i_q++) { // check if they are consistently
188 // close
189 if (verbose_)
190 ldmx_log(debug)
191 << "Comparing " << charge_check[i + i_q] << " (sample "
192 << i + i_q << ") to " << charge_check[i + 4 + i_q]
193 << " (sample " << i + 4 + i_q << "), ratio is "
194 << charge_check[i + i_q] / charge_check[i + 4 + i_q];
195 // we can be generous in these crietira since we will require an
196 // unbroken suite of 8 matches to call it oscillation
197 if (fabs(charge_check[i + i_q] / charge_check[i + 4 + i_q] - 1) <
198 0.5 || // need this tolerance to be kind of large, most
199 // actual peaks won't pass it by far anyway.
200 (charge_check[i + 4 + i_q] < 0.01 &&
201 fabs(charge_check[i + i_q] / charge_check[i + 4 + i_q]) <
202 5)) // for very small numbers, one ADC difference can be a
203 // factor 3 so add some margin
204 last_match_sample = i + i_q;
205 else {
206 if (verbose_)
207 ldmx_log(debug)
208 << "Oscillation check for channel " << digi.getChanID()
209 << " breaking at time sample " << i + i_q;
210 do_break = true; // break outer loop
211 break; // break this loop
212 }
213 }
214 if (do_break) {
215 break;
216 }
217 if (verbose_)
218 ldmx_log(debug) << "Current lastMatchSample " << last_match_sample;
219 if (last_match_sample - max_id >= 2 * 4) {
220 // we had at least a couple of oscillations (2nd period
221 // was fully matched by third)
222 flag_oscillation = 1; // there is another check possible later too,
223 // commented for now
224 break; // we've seen what we need to see
225 }
226 i += 3;
227 }
228 } // if positive maxQ
229 }
230 // //use the top and bottom ends of the sorted q as another oscillation
231 // catcher: we don't expect that the top values will be high and basically
232 // identical unless they are from an oscillation
233 int n_high = 0;
234 int quart_length = (int)charge.size() / 4;
235 for (int i = 4 * quart_length - 2; i >= 3 * quart_length; i--) {
236 // maxQ is already at last index
237 if (charge[i] / max_q > 0.66) n_high++;
238 }
239
240 /* used to be: >0.9, > 0.25. But this really killed MC pulses which
241 all end up in 0.90-0.94, and somewhere >0.05 (at least >0.15)
242 */
243 uint flag_spike = (max_q / out_event.getTotQ() > 0.95) ||
244 (charge[charge.size() - 2] / max_q < 0.05);
245 // skip "unnaturally" narrow hits (all charge in
246 // one sample or huge drop to second highest)
247 uint flag_plateau = (ped > 15 || n_high >= 5);
248 //( fabs(ped) > 15 ); //threshold // //skip
249 // events that have strange plateaus
250 uint flag_long_pulse = 0;
251 // easier to deal with in hit reconstruction directly. copy channel
252 // flags to hit flags and add this one there
253 uint flag_noise = (out_event.getNoise() > 3.5 || out_event.getNoise() == 0);
254 // =0 is typically from funky events but
255 // could be too harsh maybe
256 /* //let this wait for now
257 //if we have many high counts, a small event pedestal (where they weren't
258 included), and an avgQ ~ maxQ/nHigh, then this is an oscillation if (
259 (quartLength-nHigh<2 && ped<10) || fabs( maxQ/(avgQ*nHigh)-1 ) <0.1 )
260 flagOscillation=1;
261*/
262
263 /* //this seems to not catch what we want it to
264if ( fabs(chan.getPedestal()) < 15 //threshold // //skip events that have
265strange plateaus
266 // && (chan.getAvgQ()/chan.getPedestal()<0.8) //skip
267events that have strange oscillations
268 && 1 < nSampAboveThr && nSampAboveThr < 10 ) // skip one-time sample
269flips and long weird pulses
270 */
271 uint flag = flag_spike + 2 * flag_plateau + 4 * flag_long_pulse +
272 8 * flag_oscillation + 16 * flag_noise;
273 ldmx_log(debug)
274 << "Got quality flag " << flag
275 << " made up of (spike/plateau/long pulse/oscillation/noise) "
276 << flag_spike << "+" << flag_plateau << "+" << flag_long_pulse << "+"
277 << flag_oscillation << "+" << flag_noise;
278 out_event.setQualityFlag(flag);
279
280 // if (ped > 15 )
281 // continue;
282 ldmx_log(debug) << "In event " << event.getEventHeader().getEventNumber()
283 << ", set pedestal = " << out_event.getPedestal()
284 << // ped <<
285 " fC, noise = " << out_event.getNoise() << " fC for channel "
286 << out_event.getChanID();
287
288 channel_readout_events.push_back(out_event);
289 }
290 // Create the container to hold the
291 // digitized trigger scintillator hits.
292
293 event.add(output_collection_, channel_readout_events);
294 ldmx_log(debug) << "\n";
295}
296} // namespace trigscint
297
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Class that builds linearized full event readout.
Implements an event buffer system for storing event data.
Definition Event.h:42
Class which represents the process under execution.
Definition Process.h:36
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
Linearizes ADC info to charge, calculates channel pedestal and noise levels (in charge)
This class represents the linearised QIE output from the trigger scintillator, in charge (fC).
void setMedQ(const float medQ)
Set channel (linearized, charge-equiv) median charge.
void setFiberNb(const int fiberNb)
Set channel readout fiber number.
void setTimeOffset(const int timeOffset)
Set channel readout itme offset (in units of samples)
float getPedestal() const
Get the pedestal.
void setMinQ(const float minQ)
Set channel (linearized, charge-equiv) minimum charge.
void setQError(const std::vector< float > qErr)
Store charge quantization errors of all time samples.
void setNoise(const float noise)
Set channel (linearized, charge-equiv) noise.
float getNoise() const
Get the channel noise.
void setTotQ(const float totQ)
Set channel (linearized, charge-equiv) average charge.
float getTotQ() const
Get the channel totQ.
void setQualityFlag(const uint flag)
Set channel data quality flag.
void setPedestal(const float pedestal)
Set channel (linearized.
void setAvgQ(const float avgQ)
Set channel (linearized, charge-equiv) average charge.
void setMaxQ(const float maxQ)
Set channel (linearized, charge-equiv) maximum charge.
int getFiberNb() const
Get the channel fiberNb.
void setQ(const std::vector< float > q)
Store charges of all time samples.
void setEarlyPedestal(const float earlyPed)
Set channel (linearized.
class for simulating QIE chip output
Definition SimQIE.h:17
float adc2Q(int ADC)
Converting ADC back to charge.
Definition SimQIE.cxx:52
float qErr(float Q)
Quantization error.
Definition SimQIE.cxx:71
class for storing QIE output
void setTDC(const std::vector< int > tdc)
Store tdcs of all time samples.
void setChanID(const int chanid)
Store the channel ID.
void setTimeSinceSpill(const uint32_t timeSpill)
Store the event time since spill counter.
int getChanID() const
Get channel ID.
void setElecID(const int elecid)
Store the electronics ID.
void setADC(const std::vector< int > adc)
Store adcs of all time samples.