Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Sensor.cc
Go to the documentation of this file.
1#include <iostream>
2#include <fstream>
3#include <cmath>
4
5#include "Sensor.hh"
7#include "Plotting.hh"
8#include "Numerics.hh"
10
11namespace Garfield {
12
13double Sensor::m_signalConversion = ElementaryCharge;
14
16 : m_lastComponent(-1),
17 m_nTimeBins(200),
18 m_tStart(0.),
19 m_tStep(10.),
20 m_nEvents(0),
21 m_hasTransferFunction(false),
22 m_fTransfer(0),
23 m_hasNoiseFunction(false),
24 m_fNoise(0),
25 m_hasUserArea(false),
26 m_xMinUser(0.),
27 m_yMinUser(0.),
28 m_zMinUser(0.),
29 m_xMaxUser(0.),
30 m_yMaxUser(0.),
31 m_zMaxUser(0.),
32 m_debug(false) {
33
34 m_className = "Sensor";
35
36}
37
38ComponentBase* Sensor::GetComponent(const unsigned int i) {
39
40 if (i >= m_components.size()) {
41 std::cerr << m_className << "::GetComponent:\n";
42 std::cerr << " Component " << i << " does not exist.\n";
43 return NULL;
44 };
45 return m_components[i].comp;
46}
47
48void Sensor::ElectricField(const double x, const double y, const double z,
49 double& ex, double& ey, double& ez, double& v,
50 Medium*& medium, int& status) {
51
52 ex = ey = ez = v = 0.;
53 status = -10;
54 medium = NULL;
55 double fx, fy, fz, p;
56 Medium* med = NULL;
57 int stat;
58 // Add up electric field contributions from all components.
59 const unsigned int nComponents = m_components.size();
60 for (unsigned int i = 0; i < nComponents; ++i) {
61 m_components[i].comp->ElectricField(x, y, z, fx, fy, fz, p, med, stat);
62 if (status != 0) {
63 status = stat;
64 medium = med;
65 }
66 if (stat == 0) {
67 ex += fx;
68 ey += fy;
69 ez += fz;
70 v += p;
71 }
72 }
73}
74
75void Sensor::ElectricField(const double x, const double y, const double z,
76 double& ex, double& ey, double& ez, Medium*& medium,
77 int& status) {
78
79 ex = ey = ez = 0.;
80 status = -10;
81 medium = NULL;
82 double fx, fy, fz;
83 Medium* med = NULL;
84 int stat;
85 // Add up electric field contributions from all components.
86 const unsigned int nComponents = m_components.size();
87 for (unsigned int i = 0; i < nComponents; ++i) {
88 m_components[i].comp->ElectricField(x, y, z, fx, fy, fz, med, stat);
89 if (status != 0) {
90 status = stat;
91 medium = med;
92 }
93 if (stat == 0) {
94 ex += fx;
95 ey += fy;
96 ez += fz;
97 }
98 }
99}
100
101void Sensor::MagneticField(const double x, const double y, const double z,
102 double& bx, double& by, double& bz, int& status) {
103
104 bx = by = bz = 0.;
105 double fx, fy, fz;
106 // Add up contributions.
107 const unsigned int nComponents = m_components.size();
108 for (unsigned int i = 0; i < nComponents; ++i) {
109 m_components[i].comp->MagneticField(x, y, z, fx, fy, fz, status);
110 if (status != 0) continue;
111 bx += fx;
112 by += fy;
113 bz += fz;
114 }
115}
116
117void Sensor::WeightingField(const double x, const double y, const double z,
118 double& wx, double& wy, double& wz,
119 const std::string& label) {
120
121 wx = wy = wz = 0.;
122 double fx = 0., fy = 0., fz = 0.;
123 // Add up field contributions from all components.
124 const unsigned int nElectrodes = m_electrodes.size();
125 for (unsigned int i = 0; i < nElectrodes; ++i) {
126 if (m_electrodes[i].label == label) {
127 fx = fy = fz = 0.;
128 m_electrodes[i].comp->WeightingField(x, y, z, fx, fy, fz, label);
129 wx += fx;
130 wy += fy;
131 wz += fz;
132 }
133 }
134}
135
136double Sensor::WeightingPotential(const double x, const double y,
137 const double z, const std::string& label) {
138
139 double v = 0.;
140 // Add up contributions from all components.
141 const unsigned int nElectrodes = m_electrodes.size();
142 for (unsigned int i = 0; i < nElectrodes; ++i) {
143 if (m_electrodes[i].label == label) {
144 v += m_electrodes[i].comp->WeightingPotential(x, y, z, label);
145 }
146 }
147 return v;
148}
149
150bool Sensor::GetMedium(const double x, const double y, const double z,
151 Medium*& m) {
152
153 m = NULL;
154
155 // Make sure there is at least one component.
156 if (m_lastComponent < 0) return false;
157
158 // Check if we are still in the same component as in the previous call.
159 m = m_components[m_lastComponent].comp->GetMedium(x, y, z);
160 // Cross-check that the medium is defined.
161 if (m) return true;
162
163 const unsigned int nComponents = m_components.size();
164 for (unsigned int i = 0; i < nComponents; ++i) {
165 m = m_components[i].comp->GetMedium(x, y, z);
166 // Cross-check that the medium is defined.
167 if (m) {
168 m_lastComponent = i;
169 return true;
170 }
171 }
172 return false;
173}
174
176
177 if (!GetBoundingBox(m_xMinUser, m_yMinUser, m_zMinUser, m_xMaxUser,
178 m_yMaxUser, m_zMaxUser)) {
179 std::cerr << m_className << "::SetArea:\n";
180 std::cerr << " Bounding box is not known.\n";
181 return false;
182 }
183
184 std::cout << m_className << "::SetArea:\n";
185 std::cout << " " << m_xMinUser << " < x [cm] < " << m_xMaxUser << "\n";
186 std::cout << " " << m_yMinUser << " < y [cm] < " << m_yMaxUser << "\n";
187 std::cout << " " << m_zMinUser << " < z [cm] < " << m_zMaxUser << "\n";
188 if (std::isinf(m_xMinUser) || std::isinf(m_xMaxUser)) {
189 std::cerr << m_className << "::SetArea:\n";
190 std::cerr << " Warning: infinite x-range\n";
191 }
192 if (std::isinf(m_yMinUser) || std::isinf(m_yMaxUser)) {
193 std::cerr << m_className << "::SetArea:\n";
194 std::cerr << " Warning: infinite x-range\n";
195 }
196 if (std::isinf(m_zMinUser) || std::isinf(m_zMaxUser)) {
197 std::cerr << m_className << "::SetArea:\n";
198 std::cerr << " Warning: infinite x-range\n";
199 }
200 m_hasUserArea = true;
201 return true;
202}
203
204bool Sensor::SetArea(const double xmin, const double ymin, const double zmin,
205 const double xmax, const double ymax, const double zmax) {
206
207 if (fabs(xmax - xmin) < Small || fabs(ymax - ymin) < Small ||
208 fabs(zmax - zmin) < Small) {
209 std::cerr << m_className << "::SetArea:\n";
210 std::cerr << " Invalid range.\n";
211 return false;
212 }
213
214 m_xMinUser = xmin;
215 m_yMinUser = ymin;
216 m_zMinUser = zmin;
217 m_xMaxUser = xmax;
218 m_yMaxUser = ymax;
219 m_zMaxUser = zmax;
220
221 if (xmin > xmax) {
222 m_xMinUser = xmax;
223 m_xMaxUser = xmin;
224 }
225 if (ymin > ymax) {
226 m_yMinUser = ymax;
227 m_yMaxUser = ymin;
228 }
229 if (zmin > zmax) {
230 m_zMinUser = zmax;
231 m_zMaxUser = zmin;
232 }
233 m_hasUserArea = true;
234 return true;
235}
236
237bool Sensor::GetArea(double& xmin, double& ymin, double& zmin, double& xmax,
238 double& ymax, double& zmax) {
239
240 if (m_hasUserArea) {
241 xmin = m_xMinUser;
242 ymin = m_yMinUser;
243 zmin = m_zMinUser;
244 xmax = m_xMaxUser;
245 ymax = m_yMaxUser;
246 zmax = m_zMaxUser;
247 return true;
248 }
249
250 // User area bounds are not (yet) defined.
251 // Get the bounding box of the sensor.
252 if (!SetArea()) return false;
253
254 xmin = m_xMinUser;
255 ymin = m_yMinUser;
256 zmin = m_zMinUser;
257 xmax = m_xMaxUser;
258 ymax = m_yMaxUser;
259 zmax = m_zMaxUser;
260
261 return true;
262}
263
264bool Sensor::IsInArea(const double x, const double y, const double z) {
265
266 if (!m_hasUserArea) {
267 if (!SetArea()) {
268 std::cerr << m_className << "::IsInArea:\n";
269 std::cerr << " User area cannot be established.\n";
270 return false;
271 }
272 m_hasUserArea = true;
273 }
274
275 if (x >= m_xMinUser && x <= m_xMaxUser && y >= m_yMinUser &&
276 y <= m_yMaxUser && z >= m_zMinUser && z <= m_zMaxUser) {
277 return true;
278 }
279
280 if (m_debug) {
281 std::cout << m_className << "::IsInArea:\n" << std::endl;
282 std::cout << " (" << x << ", " << y << ", " << z << ") "
283 << " is outside.\n";
284 }
285 return false;
286}
287
288bool Sensor::IsWireCrossed(const double x0, const double y0, const double z0,
289 const double x1, const double y1, const double z1,
290 double& xc, double& yc, double& zc) {
291 const unsigned int nComponents = m_components.size();
292 for (unsigned int i = 0; i < nComponents; ++i) {
293 if (m_components[i].comp->IsWireCrossed(x0, y0, z0,
294 x1, y1, z1, xc, yc, zc)) {
295 return true;
296 }
297 }
298 return false;
299}
300
301bool Sensor::IsInTrapRadius(const double q0, const double x0,
302 const double y0, double z0, double& xw,
303 double& yw, double& rw) {
304
305 const unsigned int nComponents = m_components.size();
306 for (unsigned int i = 0; i < nComponents; ++i) {
307 if (m_components[i].comp->IsInTrapRadius(q0, x0, y0, z0, xw, yw, rw)) {
308 return true;
309 }
310 }
311 return false;
312}
313
315
316 if (!comp) {
317 std::cerr << m_className << "::AddComponent:\n Null pointer.\n";
318 return;
319 }
320
321 component newComponent;
322 newComponent.comp = comp;
323 m_components.push_back(newComponent);
324 if (m_components.size() == 1) m_lastComponent = 0;
325}
326
327void Sensor::AddElectrode(ComponentBase* comp, const std::string& label) {
328
329 if (!comp) {
330 std::cerr << m_className << "::AddElectrode:\n Null pointer.\n";
331 return;
332 }
333 const unsigned int nElectrodes = m_electrodes.size();
334 for (unsigned int i = 0; i < nElectrodes; ++i) {
335 if (m_electrodes[i].label == label) {
336 std::cout << m_className << "::AddElectrode:\n";
337 std::cout << " Warning: An electrode with label \"" << label
338 << "\" exists already.\n";
339 std::cout << " Weighting fields will be summed up.\n";
340 break;
341 }
342 }
343
344 electrode newElectrode;
345 newElectrode.comp = comp;
346 newElectrode.label = label;
347 m_electrodes.push_back(newElectrode);
348 m_electrodes.back().signal.resize(m_nTimeBins);
349 m_electrodes.back().electronsignal.resize(m_nTimeBins);
350 m_electrodes.back().ionsignal.resize(m_nTimeBins);
351 std::cout << m_className << "::AddElectrode:\n";
352 std::cout << " Added readout electrode \"" << label << "\".\n";
353 std::cout << " All signals are reset.\n";
354 ClearSignal();
355}
356
358
359 m_components.clear();
360 m_lastComponent = -1;
361 m_electrodes.clear();
362 m_nTimeBins = 200;
363 m_tStart = 0.;
364 m_tStep = 10.;
365 m_nEvents = 0;
366 m_hasUserArea = false;
367}
368
369bool Sensor::GetVoltageRange(double& vmin, double& vmax) {
370
371 // We don't know the range yet.
372 bool set = false;
373 // Loop over the components.
374 double umin, umax;
375 const unsigned int nComponents = m_components.size();
376 for (unsigned int i = 0; i < nComponents; ++i) {
377 if (!m_components[i].comp->GetVoltageRange(umin, umax)) continue;
378 if (set) {
379 if (umin < vmin) vmin = umin;
380 if (umax > vmax) vmax = umax;
381 } else {
382 vmin = umin;
383 vmax = umax;
384 set = true;
385 }
386 }
387
388 // Warn if we still don't know the range.
389 if (!set) {
390 std::cerr << m_className << "::GetVoltageRange:\n";
391 std::cerr << " Sensor voltage range not known.\n";
392 vmin = vmax = 0.;
393 return false;
394 }
395
396 if (m_debug) {
397 std::cout << m_className << "::GetVoltageRange:\n";
398 std::cout << " Voltage range " << vmin << " < V < " << vmax << ".\n";
399 }
400 return true;
401}
402
404
405 const unsigned int nElectrodes = m_electrodes.size();
406 for (unsigned int i = 0; i < nElectrodes; ++i) {
407 m_electrodes[i].charge = 0.;
408 for (int j = m_nTimeBins; j--;) {
409 m_electrodes[i].signal[j] = 0.;
410 m_electrodes[i].electronsignal[j] = 0.;
411 m_electrodes[i].ionsignal[j] = 0.;
412 }
413 }
414 m_nEvents = 0;
415}
416
417void Sensor::AddSignal(const double q, const double t, const double dt,
418 const double x, const double y, const double z,
419 const double vx, const double vy, const double vz) {
420 // Get the time bin.
421 if (t < m_tStart || dt <= 0.) {
422 if (m_debug) {
423 std::cerr << m_className << "::AddSignal:\n";
424 if (t < m_tStart) std::cerr << " Time " << t << " out of range.\n";
425 if (dt <= 0.) std::cerr << " Time step < 0.\n";
426 }
427 return;
428 }
429 const int bin = int((t - m_tStart) / m_tStep);
430 // Check if the starting time is outside the range
431 if (bin < 0 || bin >= (int)m_nTimeBins) {
432 if (m_debug) {
433 std::cerr << m_className << "::AddSignal:\n";
434 std::cerr << " Bin " << bin << " out of range.\n";
435 }
436 return;
437 }
438 if (m_nEvents <= 0) m_nEvents = 1;
439
440 double wx = 0., wy = 0., wz = 0.;
441 if (m_debug) {
442 std::cout << m_className << "::AddSignal:\n";
443 std::cout << " Time: " << t << "\n";
444 std::cout << " Step: " << dt << "\n";
445 std::cout << " Charge: " << q << "\n";
446 std::cout << " Velocity: (" << vx << ", " << vy << ", " << vz << ")\n";
447 }
448 const unsigned int nElectrodes = m_electrodes.size();
449 for (unsigned int i = 0; i < nElectrodes; ++i) {
450 // Calculate the weighting field for this electrode
451 m_electrodes[i].comp->WeightingField(x, y, z, wx, wy, wz,
452 m_electrodes[i].label);
453 // Calculate the induced current
454 const double cur = -q * (wx * vx + wy * vy + wz * vz);
455 if (m_debug) {
456 std::cout << " Electrode " << m_electrodes[i].label << ":\n";
457 std::cout << " Weighting field: (" << wx << ", " << wy << ", " << wz
458 << ")\n";
459 std::cout << " Induced charge: " << cur* dt << "\n";
460 }
461 double delta = m_tStart + (bin + 1) * m_tStep - t;
462 // Check if the provided timestep extends over more than one time bin
463 if (dt > delta) {
464 m_electrodes[i].signal[bin] += cur * delta;
465 if (q < 0) {
466 m_electrodes[i].electronsignal[bin] += cur * delta;
467 } else {
468 m_electrodes[i].ionsignal[bin] += cur * delta;
469 }
470 delta = dt - delta;
471 unsigned int j = 1;
472 while (delta > m_tStep && bin + j < m_nTimeBins) {
473 m_electrodes[i].signal[bin + j] += cur * m_tStep;
474 if (q < 0) {
475 m_electrodes[i].electronsignal[bin + j] += cur * m_tStep;
476 } else {
477 m_electrodes[i].ionsignal[bin + j] += cur * m_tStep;
478 }
479 delta -= m_tStep;
480 ++j;
481 }
482 if (bin + j < m_nTimeBins) {
483 m_electrodes[i].signal[bin + j] += cur * delta;
484 if (q < 0) {
485 m_electrodes[i].electronsignal[bin + j] += cur * delta;
486 } else {
487 m_electrodes[i].ionsignal[bin + j] += cur * delta;
488 }
489 }
490 } else {
491 m_electrodes[i].signal[bin] += cur * dt;
492 if (q < 0) {
493 m_electrodes[i].electronsignal[bin] += cur * dt;
494 } else {
495 m_electrodes[i].ionsignal[bin] += cur * dt;
496 }
497 }
498 }
499}
500
501void Sensor::AddInducedCharge(const double q, const double x0, const double y0,
502 const double z0, const double x1, const double y1,
503 const double z1) {
504
505 if (m_debug) std::cout << m_className << "::AddInducedCharge:\n";
506 double w0 = 0., w1 = 0.;
507 const unsigned int nElectrodes = m_electrodes.size();
508 for (unsigned int i = 0; i < nElectrodes; ++i) {
509 // Calculate the weighting potential at the starting point.
510 w0 = m_electrodes[i]
511 .comp->WeightingPotential(x0, y0, z0, m_electrodes[i].label);
512 // Calculate the weighting potential at the end point.
513 w1 = m_electrodes[i]
514 .comp->WeightingPotential(x1, y1, z1, m_electrodes[i].label);
515 m_electrodes[i].charge += q * (w1 - w0);
516 if (m_debug) {
517 std::cout << " Electrode " << m_electrodes[i].label << ":\n";
518 std::cout << " Weighting potential at (" << x0 << ", " << y0 << ", "
519 << z0 << "): " << w0 << "\n";
520 std::cout << " Weighting potential at (" << x1 << ", " << y1 << ", "
521 << z1 << "): " << w1 << "\n";
522 std::cout << " Induced charge: " << m_electrodes[i].charge << "\n";
523 }
524 }
525}
526
527void Sensor::SetTimeWindow(const double tstart, const double tstep,
528 const unsigned int nsteps) {
529
530 m_tStart = tstart;
531 if (tstep <= 0.) {
532 std::cerr << m_className << "::SetTimeWindow:\n";
533 std::cerr << " Starting time out of range.\n";
534 } else {
535 m_tStep = tstep;
536 }
537
538 if (nsteps == 0) {
539 std::cerr << m_className << "::SetTimeWindow:\n";
540 std::cerr << " Number of time bins out of range.\n";
541 } else {
542 m_nTimeBins = nsteps;
543 }
544
545 if (m_debug) {
546 std::cout << m_className << "::SetTimeWindow:\n";
547 std::cout << " " << m_tStart << " < t [ns] < "
548 << m_tStart + m_nTimeBins* m_tStep << "\n";
549 std::cout << " Step size: " << m_tStep << " ns\n";
550 }
551
552 std::cout << m_className << "::SetTimeWindow:\n";
553 std::cout << " Resetting all signals.\n";
554 const unsigned int nElectrodes = m_electrodes.size();
555 for (unsigned int i = 0; i < nElectrodes; ++i) {
556 m_electrodes[i].signal.assign(m_nTimeBins, 0.);
557 m_electrodes[i].electronsignal.assign(m_nTimeBins, 0.);
558 m_electrodes[i].ionsignal.assign(m_nTimeBins, 0.);
559 }
560 m_nEvents = 0;
561}
562
563double Sensor::GetElectronSignal(const std::string& label,
564 const unsigned int bin) {
565
566 if (m_nEvents == 0) return 0.;
567 if (bin >= m_nTimeBins) return 0.;
568 double sig = 0.;
569 const unsigned int nElectrodes = m_electrodes.size();
570 for (unsigned int i = 0; i < nElectrodes; ++i) {
571 if (m_electrodes[i].label == label)
572 sig += m_electrodes[i].electronsignal[bin];
573 }
574 if (m_debug) {
575 std::cout << m_className << "::GetElectronSignal:\n";
576 std::cout << " Electrode: " << label << "\n";
577 std::cout << " Bin: " << bin << "\n";
578 std::cout << " ElectronSignal: " << sig / m_tStep << "\n";
579 }
580 return m_signalConversion * sig / (m_nEvents * m_tStep);
581}
582
583double Sensor::GetIonSignal(const std::string& label, const unsigned int bin) {
584
585 if (m_nEvents == 0) return 0.;
586 if (bin >= m_nTimeBins) return 0.;
587 double sig = 0.;
588 const unsigned int nElectrodes = m_electrodes.size();
589 for (unsigned int i = 0; i < nElectrodes; ++i) {
590 if (m_electrodes[i].label == label) sig += m_electrodes[i].ionsignal[bin];
591 }
592 if (m_debug) {
593 std::cout << m_className << "::GetIonSignal:\n";
594 std::cout << " Electrode: " << label << "\n";
595 std::cout << " Bin: " << bin << "\n";
596 std::cout << " IonSignal: " << sig / m_tStep << "\n";
597 }
598 return m_signalConversion * sig / (m_nEvents * m_tStep);
599}
600
601double Sensor::GetSignal(const std::string& label, const unsigned int bin) {
602
603 if (m_nEvents == 0) return 0.;
604 if (bin >= m_nTimeBins) return 0.;
605 double sig = 0.;
606 const unsigned int nElectrodes = m_electrodes.size();
607 for (unsigned int i = 0; i < nElectrodes; ++i) {
608 if (m_electrodes[i].label == label) sig += m_electrodes[i].signal[bin];
609 }
610 if (m_debug) {
611 std::cout << m_className << "::GetSignal:\n";
612 std::cout << " Electrode: " << label << "\n";
613 std::cout << " Bin: " << bin << "\n";
614 std::cout << " Signal: " << sig / m_tStep << "\n";
615 }
616 return m_signalConversion * sig / (m_nEvents * m_tStep);
617}
618
619double Sensor::GetInducedCharge(const std::string& label) {
620
621 if (m_nEvents == 0) return 0.;
622 double charge = 0.;
623 const unsigned int nElectrodes = m_electrodes.size();
624 for (unsigned int i = 0; i < nElectrodes; ++i) {
625 if (m_electrodes[i].label == label) charge += m_electrodes[i].charge;
626 }
627 if (m_debug) {
628 std::cout << m_className << "::GetInducedCharge:\n";
629 std::cout << " Electrode: " << label << "\n";
630 std::cout << " Charge: " << charge / m_tStep << "\n";
631 }
632
633 return charge / m_nEvents;
634}
635
636void Sensor::SetTransferFunction(double (*f)(double t)) {
637
638 if (!f) {
639 std::cerr << m_className << "::SetTransferFunction:\n Null pointer.\n";
640 return;
641 }
642 m_fTransfer = f;
643 m_hasTransferFunction = true;
644 m_transferFunctionTimes.clear();
645 m_transferFunctionValues.clear();
646}
647
648void Sensor::SetTransferFunction(const std::vector<double>& times,
649 const std::vector<double>& values) {
650
651 if (times.empty() || values.empty()) {
652 std::cerr << m_className << "::SetTransferFunction:\n";
653 std::cerr << " Time and value vectors must not be empty.\n";
654 return;
655 } else if (times.size() != values.size()) {
656 std::cerr << m_className << "::SetTransferFunction:\n";
657 std::cerr << " Time and value vectors must have same size.\n";
658 return;
659 }
660 m_transferFunctionTimes = times;
661 m_transferFunctionValues = values;
662 m_fTransfer = NULL;
663 m_hasTransferFunction = true;
664}
665
666double Sensor::InterpolateTransferFunctionTable(double t) {
667
668 if (m_transferFunctionTimes.empty() || m_transferFunctionValues.empty()) {
669 return 0.;
670 }
671 // Don't extrapolate beyond the range defined in the table.
672 if (t < m_transferFunctionTimes.front() ||
673 t > m_transferFunctionTimes.back()) {
674 return 0.;
675 }
676 // Find the proper interval in the table.
677 int iLow = 0;
678 int iUp = m_transferFunctionTimes.size() - 1;
679 int iM;
680 while (iUp - iLow > 1) {
681 iM = (iUp + iLow) >> 1;
682 if (t >= m_transferFunctionTimes[iM]) {
683 iLow = iM;
684 } else {
685 iUp = iM;
686 }
687 }
688 // Linear interpolation.
689 return m_transferFunctionValues[iLow] +
690 (t - m_transferFunctionTimes[iLow]) *
691 (m_transferFunctionValues[iUp] - m_transferFunctionValues[iLow]) /
692 (m_transferFunctionTimes[iUp] - m_transferFunctionTimes[iLow]);
693}
694
695double Sensor::GetTransferFunction(const double t) {
696
697 if (!m_hasTransferFunction) return 0.;
698 if (m_fTransfer) return m_fTransfer(t);
699 return InterpolateTransferFunctionTable(t);
700}
701
703
704 if (!m_hasTransferFunction) {
705 std::cerr << m_className << "::ConvoluteSignal:\n";
706 std::cerr << " No transfer function available.\n";
707 return false;
708 }
709 if (m_nEvents == 0) {
710 std::cerr << m_className << "::ConvoluteSignal:\n";
711 std::cerr << " No signals present.\n";
712 return false;
713 }
714
715 // Set the range where the transfer function is valid.
716 double cnvMin = 0.;
717 double cnvMax = 1.e10;
718
719 std::vector<double> cnvTab(2 * m_nTimeBins, 0.);
720 int iOffset = m_nTimeBins;
721 // Evaluate the transfer function.
722 for (unsigned int i = 0; i < m_nTimeBins; ++i) {
723 // Negative time part.
724 double t = -i * m_tStep;
725 if (t < cnvMin || t > cnvMax) {
726 cnvTab[iOffset - i] = 0.;
727 } else if (m_fTransfer) {
728 cnvTab[iOffset - i] = m_fTransfer(t);
729 } else {
730 cnvTab[iOffset - i] = InterpolateTransferFunctionTable(t);
731 }
732 if (i < 1) continue;
733 // Positive time part.
734 t = i * m_tStep;
735 if (t < cnvMin || t > cnvMax) {
736 cnvTab[iOffset + i] = 0.;
737 } else if (m_fTransfer) {
738 cnvTab[iOffset + i] = m_fTransfer(t);
739 } else {
740 cnvTab[iOffset + i] = InterpolateTransferFunctionTable(t);
741 }
742 }
743
744 std::vector<double> tmpSignal(m_nTimeBins, 0.);
745 // Loop over all electrodes.
746 const unsigned int nElectrodes = m_electrodes.size();
747 for (unsigned int i = 0; i < nElectrodes; ++i) {
748 for (unsigned int j = 0; j < m_nTimeBins; ++j) {
749 tmpSignal[j] = 0.;
750 for (unsigned int k = 0; k < m_nTimeBins; ++k) {
751 tmpSignal[j] +=
752 m_tStep * cnvTab[iOffset + j - k] * m_electrodes[i].signal[k];
753 }
754 }
755 m_electrodes[i].signal = tmpSignal;
756 }
757 return true;
758}
759
761
762 if (m_nEvents == 0) {
763 std::cerr << m_className << "::IntegrateSignal:\n";
764 std::cerr << " No signals present.\n";
765 return false;
766 }
767
768 const unsigned int nElectrodes = m_electrodes.size();
769 for (unsigned int i = 0; i < nElectrodes; ++i) {
770 for (unsigned int j = 0; j < m_nTimeBins; ++j) {
771 m_electrodes[i].signal[j] *= m_tStep;
772 m_electrodes[i].electronsignal[j] *= m_tStep;
773 m_electrodes[i].ionsignal[j] *= m_tStep;
774 if (j > 0) {
775 m_electrodes[i].signal[j] += m_electrodes[i].signal[j - 1];
776 m_electrodes[i].electronsignal[j] +=
777 m_electrodes[i].electronsignal[j - 1];
778 m_electrodes[i].ionsignal[j] += m_electrodes[i].ionsignal[j - 1];
779 }
780 }
781 }
782 return true;
783}
784
785void Sensor::SetNoiseFunction(double (*f)(double t)) {
786
787 if (f == 0) {
788 std::cerr << m_className << "::SetNoiseFunction:\n";
789 std::cerr << " Function pointer is null.\n";
790 return;
791 }
792 m_fNoise = f;
793 m_hasNoiseFunction = true;
794}
795
796void Sensor::AddNoise(const bool total, const bool electron, const bool ion) {
797
798 if (!m_hasNoiseFunction) {
799 std::cerr << m_className << "::AddNoise:\n";
800 std::cerr << " Noise function is not defined.\n";
801 return;
802 }
803 if (m_nEvents == 0) m_nEvents = 1;
804
805 const unsigned int nElectrodes = m_electrodes.size();
806 for (unsigned int i = 0; i < nElectrodes; ++i) {
807 for (unsigned int j = 0; j < m_nTimeBins; ++j) {
808 const double t = m_tStart + (j + 0.5) * m_tStep;
809 const double noise = m_fNoise(t);
810 if (total) m_electrodes[i].signal[j] += noise;
811 if (electron) m_electrodes[i].electronsignal[j] += noise;
812 if (ion) m_electrodes[i].ionsignal[j] += noise;
813 }
814 }
815}
816
818 const std::string& label, int& n) {
819
820 // Reset the list of threshold crossings.
821 m_thresholdCrossings.clear();
822 m_thresholdLevel = thr;
823
824 // Set the interpolation order.
825 int iOrder = 1;
826
827 if (m_nEvents == 0) {
828 std::cerr << m_className << "::ComputeThresholdCrossings:\n";
829 std::cerr << " No signals present.\n";
830 return false;
831 }
832
833 // Compute the total signal.
834 std::vector<double> signal(m_nTimeBins, 0.);
835 // Loop over the electrodes.
836 bool foundLabel = false;
837 const unsigned int nElectrodes = m_electrodes.size();
838 for (unsigned int j = 0; j < nElectrodes; ++j) {
839 if (m_electrodes[j].label == label) {
840 foundLabel = true;
841 for (unsigned int i = 0; i < m_nTimeBins; ++i) {
842 signal[i] += m_electrodes[j].signal[i];
843 }
844 }
845 }
846 if (!foundLabel) {
847 std::cerr << m_className << "::ComputeThresholdCrossings:\n";
848 std::cerr << " Electrode " << label << " not found.\n";
849 return false;
850 }
851 const double scale = m_signalConversion / (m_nEvents * m_tStep);
852 for (unsigned int i = 0; i < m_nTimeBins; ++i) signal[i] *= scale;
853
854 // Establish the range.
855 double vMin = signal[0];
856 double vMax = signal[0];
857 for (unsigned int i = 0; i < m_nTimeBins; ++i) {
858 if (signal[i] < vMin) vMin = signal[i];
859 if (signal[i] > vMax) vMax = signal[i];
860 }
861 if (thr < vMin && thr > vMax) {
862 if (m_debug) {
863 std::cout << m_className << "::ComputeThresholdCrossings:\n";
864 std::cout << " Threshold outside the range [" << vMin << ", " << vMax
865 << "]\n";
866 }
867 return true;
868 }
869
870 // Check for rising edges.
871 bool rise = true;
872 bool fall = false;
873
874 while (rise || fall) {
875 if (m_debug) {
876 if (rise) {
877 std::cout << m_className << "::ComputeThresholdCrossings:\n";
878 std::cout << " Hunting for rising edges.\n";
879 } else if (fall) {
880 std::cout << m_className << "::ComputeThresholdCrossings:\n";
881 std::cout << " Hunting for falling edges.\n";
882 }
883 }
884 // Initialise the vectors.
885 std::vector<double> times;
886 std::vector<double> values;
887 times.push_back(m_tStart);
888 values.push_back(signal[0]);
889 int nValues = 1;
890 // Scan the signal.
891 for (unsigned int i = 1; i < m_nTimeBins; ++i) {
892 // Compute the vector element.
893 const double tNew = m_tStart + i * m_tStep;
894 const double vNew = signal[i];
895 // If still increasing or decreasing, add to the vector.
896 if ((rise && vNew > values.back()) || (fall && vNew < values.back())) {
897 times.push_back(tNew);
898 values.push_back(vNew);
899 ++nValues;
900 // Otherwise see whether we crossed the threshold level.
901 } else if ((values[0] - thr) * (thr - values.back()) >= 0. &&
902 nValues > 1 && ((rise && values.back() > values[0]) ||
903 (fall && values.back() < values[0]))) {
904 // Compute the crossing time.
905 double tcr = Numerics::Divdif(times, values, nValues, thr, iOrder);
906 thresholdCrossing newCrossing;
907 newCrossing.time = tcr;
908 newCrossing.rise = rise;
909 m_thresholdCrossings.push_back(newCrossing);
910 times.clear();
911 values.clear();
912 times.push_back(tNew);
913 values.push_back(vNew);
914 nValues = 1;
915 } else {
916 // No crossing, simply reset the vector.
917 times.clear();
918 values.clear();
919 times.push_back(tNew);
920 values.push_back(vNew);
921 nValues = 1;
922 }
923 }
924 // Check the final vector.
925 if ((values[0] - thr) * (thr - values.back()) >= 0. && nValues > 1 &&
926 ((rise && values.back() > values[0]) ||
927 (fall && values.back() < values[0]))) {
928 double tcr = Numerics::Divdif(times, values, nValues, thr, iOrder);
929 thresholdCrossing newCrossing;
930 newCrossing.time = tcr;
931 newCrossing.rise = rise;
932 m_thresholdCrossings.push_back(newCrossing);
933 }
934 if (rise) {
935 rise = false;
936 fall = true;
937 } else if (fall) {
938 rise = fall = false;
939 }
940 }
941 n = m_thresholdCrossings.size();
942
943 if (m_debug) {
944 std::cout << m_className << "::ComputeThresholdCrossings:\n";
945 std::cout << " Found " << n << " crossings.\n";
946 if (n > 0) std::cout << " Time [ns] Direction\n";
947 for (int i = 0; i < n; ++i) {
948 std::cout << " " << m_thresholdCrossings[i].time << " ";
949 if (m_thresholdCrossings[i].rise) {
950 std::cout << "rising\n";
951 } else {
952 std::cout << "falling\n";
953 }
954 }
955 }
956 // Seems to have worked.
957 return true;
958}
959
960bool Sensor::GetThresholdCrossing(const unsigned int i, double& time,
961 double& level, bool& rise) const {
962
963 level = m_thresholdLevel;
964
965 if (i >= m_thresholdCrossings.size()) {
966 std::cerr << m_className << "::GetThresholdCrossing:\n";
967 std::cerr << " Index (" << i << ") out of range.\n";
968 time = m_tStart + m_nTimeBins * m_tStep;
969 return false;
970 }
971
972 time = m_thresholdCrossings[i].time;
973 rise = m_thresholdCrossings[i].rise;
974 return true;
975}
976
977bool Sensor::GetBoundingBox(double& xmin, double& ymin, double& zmin,
978 double& xmax, double& ymax, double& zmax) {
979
980 // We don't know the range yet
981 bool set = false;
982 // Loop over the fields
983 double x0, y0, z0, x1, y1, z1;
984 const unsigned int nComponents = m_components.size();
985 for (unsigned int i = 0; i < nComponents; ++i) {
986 if (!m_components[i].comp->GetBoundingBox(x0, y0, z0, x1, y1, z1)) continue;
987 if (set) {
988 if (x0 < xmin) xmin = x0;
989 if (y0 < ymin) ymin = y0;
990 if (z0 < zmin) zmin = z0;
991 if (x1 > xmax) xmax = x1;
992 if (y1 > ymax) ymax = y1;
993 if (z1 > zmax) zmax = z1;
994 } else {
995 xmin = x0;
996 ymin = y0;
997 zmin = z0;
998 xmax = x1;
999 ymax = y1;
1000 zmax = z1;
1001 set = true;
1002 }
1003 }
1004
1005 // Warn if we still don't know the range
1006 if (!set) {
1007 std::cerr << m_className << "::GetBoundingBox:\n";
1008 std::cerr << " Sensor bounding box not known.\n";
1009 xmin = 0.;
1010 ymin = 0.;
1011 zmin = 0.;
1012 xmax = 0.;
1013 ymax = 0.;
1014 zmax = 0.;
1015 return false;
1016 }
1017
1018 if (m_debug) {
1019 std::cout << m_className << "::GetBoundingBox:\n";
1020 std::cout << " " << xmin << " < x [cm] < " << xmax << "\n";
1021 std::cout << " " << ymin << " < y [cm] < " << ymax << "\n";
1022 std::cout << " " << zmin << " < z [cm] < " << zmax << "\n";
1023 }
1024 return true;
1025}
1026}
Abstract base class for components.
Abstract base class for media.
Definition: Medium.hh:11
void MagneticField(const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
Get the magnetic field at (x, y, z).
Definition: Sensor.cc:101
double GetElectronSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:563
bool GetVoltageRange(double &vmin, double &vmax)
Return the voltage range.
Definition: Sensor.cc:369
void SetTimeWindow(const double tstart, const double tstep, const unsigned int nsteps)
Definition: Sensor.cc:527
bool IntegrateSignal()
Definition: Sensor.cc:760
bool IsInTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
Definition: Sensor.cc:301
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label)
Get the weighting field at (x, y, z).
Definition: Sensor.cc:117
double GetTransferFunction(const double t)
Definition: Sensor.cc:695
Sensor()
Constructor.
Definition: Sensor.cc:15
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&medium, int &status)
Get the drift field and potential at (x, y, z).
Definition: Sensor.cc:48
bool SetArea()
Set the user area to the default.
Definition: Sensor.cc:175
double GetSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:601
void SetNoiseFunction(double(*f)(double t))
Definition: Sensor.cc:785
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
Definition: Sensor.cc:264
void AddElectrode(ComponentBase *comp, const std::string &label)
Add an electrode.
Definition: Sensor.cc:327
double GetInducedCharge(const std::string &label)
Definition: Sensor.cc:619
void AddNoise(const bool total=true, const bool electron=false, const bool ion=false)
Definition: Sensor.cc:796
void AddSignal(const double q, const double t, const double dt, const double x, const double y, const double z, const double vx, const double vy, const double vz)
Definition: Sensor.cc:417
double WeightingPotential(const double x, const double y, const double z, const std::string &label)
Get the weighting potential at (x, y, z).
Definition: Sensor.cc:136
bool ComputeThresholdCrossings(const double thr, const std::string &label, int &n)
Definition: Sensor.cc:817
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:150
ComponentBase * GetComponent(const unsigned int componentNumber)
Definition: Sensor.cc:38
void Clear()
Remove all components, electrodes and reset the sensor.
Definition: Sensor.cc:357
void SetTransferFunction(double(*f)(double t))
Definition: Sensor.cc:636
bool GetThresholdCrossing(const unsigned int i, double &time, double &level, bool &rise) const
Definition: Sensor.cc:960
bool IsWireCrossed(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc)
Definition: Sensor.cc:288
bool ConvoluteSignal()
Definition: Sensor.cc:702
void AddComponent(ComponentBase *comp)
Add a component.
Definition: Sensor.cc:314
void ClearSignal()
Reset signals and induced charges of all electrodes.
Definition: Sensor.cc:403
void AddInducedCharge(const double q, const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Definition: Sensor.cc:501
double GetIonSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:583
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
Definition: Sensor.cc:237
double Divdif(const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
Definition: Numerics.cc:650