Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewIsochrons.cc
Go to the documentation of this file.
1#include <algorithm>
2#include <cmath>
3#include <iostream>
4#include <set>
5
6#include <TAxis.h>
7#include <TROOT.h>
8#include <TGraph.h>
9#include <TH1F.h>
10
11#include "Garfield/Sensor.hh"
12#include "Garfield/Component.hh"
14#include "Garfield/Plotting.hh"
16
17namespace {
18
19double Interpolate(const std::vector<double>& y,
20 const std::vector<double>& x, const double xx) {
21
22 const double tol = 1.e-6 * fabs(x.back() - x.front());
23 if (xx < x.front()) return y.front();
24 const auto it1 = std::upper_bound(x.cbegin(), x.cend(), xx);
25 if (it1 == x.cend()) return y.back();
26 const auto it0 = std::prev(it1);
27 const double dx = (*it1 - *it0);
28 if (dx < tol) return y[it0 - x.cbegin()];
29 const double f = (xx - *it0) / dx;
30 return y[it0 - x.cbegin()] * (1. - f) + f * y[it1 - x.cbegin()];
31}
32
33bool OnLine(const double x1, const double y1, const double x2, const double y2,
34 const double u, const double v) {
35 // Set tolerances.
36 double epsx = 1.e-10 * std::max({fabs(x1), fabs(x2), fabs(u)});
37 double epsy = 1.e-10 * std::max({fabs(y1), fabs(y2), fabs(v)});
38 epsx = std::max(1.e-10, epsx);
39 epsy = std::max(1.e-10, epsy);
40
41 if ((fabs(x1 - u) <= epsx && fabs(y1 - v) <= epsy) ||
42 (fabs(x2 - u) <= epsx && fabs(y2 - v) <= epsy)) {
43 // Point to be examined coincides with start or end.
44 return true;
45 } else if (fabs(x1 - x2) <= epsx && fabs(y1 - y2) <= epsy) {
46 // The line (x1, y1) to (x2, y2) is in fact a point.
47 return false;
48 }
49 double xc = 0., yc = 0.;
50 if (fabs(u - x1) + fabs(v - y1) < fabs(u - x2) + fabs(v - y2)) {
51 // (u, v) is nearer to (x1, y1).
52 const double dx = (x2 - x1);
53 const double dy = (y2 - y1);
54 const double xl = ((u - x1) * dx + (v - y1) * dy) / (dx * dx + dy * dy);
55 if (xl < 0.) {
56 xc = x1;
57 yc = y1;
58 } else if (xl > 1.) {
59 xc = x2;
60 yc = y2;
61 } else {
62 xc = x1 + xl * dx;
63 yc = y1 + xl * dy;
64 }
65 } else {
66 // (u, v) is nearer to (x2, y2).
67 const double dx = (x1 - x2);
68 const double dy = (y1 - y2);
69 const double xl = ((u - x2) * dx + (v - y2) * dy) / (dx * dx + dy * dy);
70 if (xl < 0.) {
71 xc = x2;
72 yc = y2;
73 } else if (xl > 1.) {
74 xc = x1;
75 yc = y1;
76 } else {
77 xc = x2 + xl * dx;
78 yc = y2 + xl * dy;
79 }
80 }
81 // See whether the point is on the line.
82 if (fabs(u - xc) < epsx && fabs(v - yc) < epsy) {
83 return true;
84 }
85 return false;
86}
87
88bool Crossing(const double x1, const double y1, const double x2,
89 const double y2, const double u1, const double v1,
90 const double u2, const double v2) {
91
92 // Check for a point of one line located on the other line.
93 if (OnLine(x1, y1, x2, y2, u1, v1) || OnLine(x1, y1, x2, y2, u2, v2) ||
94 OnLine(u1, v1, u2, v2, x1, y1) || OnLine(u1, v1, u2, v2, x2, y2)) {
95 return true;
96 }
97 // Matrix to compute the crossing point.
98 std::array<std::array<double, 2>, 2> a;
99 a[0][0] = y2 - y1;
100 a[0][1] = v2 - v1;
101 a[1][0] = x1 - x2;
102 a[1][1] = u1 - u2;
103 const double det = a[0][0] * a[1][1] - a[1][0] * a[0][1];
104 // Set tolerances.
105 double epsx = 1.e-10 * std::max({fabs(x1), fabs(x2), fabs(u1), fabs(u2)});
106 double epsy = 1.e-10 * std::max({fabs(y1), fabs(y2), fabs(v1), fabs(v2)});
107 epsx = std::max(epsx, 1.e-10);
108 epsy = std::max(epsy, 1.e-10);
109 if (fabs(det) < epsx * epsy) {
110 // Parallel, non-touching.
111 return false;
112 }
113 // Crossing, non-trivial lines: solve crossing equations.
114 const double aux = a[1][1];
115 a[1][1] = a[0][0] / det;
116 a[0][0] = aux / det;
117 a[1][0] = -a[1][0] / det;
118 a[0][1] = -a[0][1] / det;
119 // Compute crossing point.
120 const double xc = a[0][0] * (x1 * y2 - x2 * y1) + a[1][0] * (u1 * v2 - u2 * v1);
121 const double yc = a[0][1] * (x1 * y2 - x2 * y1) + a[1][1] * (u1 * v2 - u2 * v1);
122 // See whether the crossing point is on both lines.
123 if (OnLine(x1, y1, x2, y2, xc, yc) && OnLine(u1, v1, u2, v2, xc, yc)) {
124 // Intersecting lines.
125 return true;
126 }
127 // Crossing point not on both lines.
128 return false;
129}
130
131}
132
133namespace Garfield {
134
136
138 if (!s) {
139 std::cerr << m_className << "::SetSensor: Null pointer.\n";
140 return;
141 }
142
143 m_sensor = s;
144 m_component = nullptr;
145}
146
148 if (!c) {
149 std::cerr << m_className << "::SetComponent: Null pointer.\n";
150 return;
151 }
152
153 m_component = c;
154 m_sensor = nullptr;
155}
156
158
159 if (ar < 0.) {
160 std::cerr << m_className << "::SetAspectRatioSwitch: Value must be > 0.\n";
161 return;
162 }
163 m_aspectRatio = ar;
164}
165
166void ViewIsochrons::SetLoopThreshold(const double thr) {
167
168 if (thr < 0. || thr > 1.) {
169 std::cerr << m_className << "::SetLoopThreshold:\n"
170 << " Value must be between 0 and 1.\n";
171 return;
172 }
173 m_loopThreshold = thr;
174}
175
177
178 if (thr < 0. || thr > 1.) {
179 std::cerr << m_className << "::SetConnectionThreshold:\n"
180 << " Value must be between 0 and 1.\n";
181 return;
182 }
183 m_connectionThreshold = thr;
184}
185
186void ViewIsochrons::PlotIsochrons(const double tstep,
187 const std::vector<std::array<double, 3> >& points, const bool rev,
188 const bool colour, const bool markers, const bool plotDriftLines) {
189
190 if (!m_sensor && !m_component) {
191 std::cerr << m_className << "::PlotIsochrons:\n"
192 << " Neither sensor nor component are defined.\n";
193 return;
194 }
195 if (tstep <= 0.) {
196 std::cerr << m_className << "::PlotIsochrons: Time step must be > 0.\n";
197 return;
198 }
199 if (points.empty()) {
200 std::cerr << m_className << "::PlotIsochrons:\n"
201 << " No starting points provided.\n";
202 return;
203 }
204 if (!SetPlotLimits()) return;
205 auto canvas = GetCanvas();
206 canvas->cd();
207 canvas->SetTitle("Isochrons");
208 auto frame = canvas->DrawFrame(m_xMinPlot, m_yMinPlot,
210 frame->GetXaxis()->SetTitle(LabelX().c_str());
211 frame->GetYaxis()->SetTitle(LabelY().c_str());
212 canvas->Update();
213
214 //-----------------------------------------------------------------------
215 // DRFEQT - The main routine (DRFEQT) accumulates equal drift time data
216 // DRFEQP which is plotted as a set of contours in the entry DRFEQP.
217 //-----------------------------------------------------------------------
218 std::vector<std::vector<std::array<double, 3> > > driftLines;
219 std::vector<std::array<double, 3> > startPoints;
220 std::vector<std::array<double, 3> > endPoints;
221 std::vector<int> statusCodes;
222 // Accumulate drift lines.
223 ComputeDriftLines(tstep, points, driftLines, startPoints, endPoints,
224 statusCodes, rev);
225 const unsigned int nDriftLines = driftLines.size();
226 if (nDriftLines < 2) {
227 std::cerr << m_className << "::PlotIsochrons: Too few drift lines.\n";
228 return;
229 }
230 // Keep track of the largest number of contours.
231 std::size_t nContours = 0;
232 for (const auto& driftLine : driftLines) {
233 nContours = std::max(nContours, driftLine.size());
234 }
235 if (nContours == 0) {
236 std::cerr << m_className << "::PlotIsochrons: No contour lines.\n";
237 return;
238 }
239
240 std::set<int> allStats;
241 for (const auto stat : statusCodes) allStats.insert(stat);
242
243 // DRFEQP
244 if (m_debug) {
245 std::cout << m_className << "::PlotIsochrons:\n"
246 << " Drawing " << nContours << " contours, "
247 << nDriftLines << " drift lines.\n";
248 std::printf(" Connection threshold: %10.3f\n", m_connectionThreshold);
249 std::printf(" Aspect ratio threshold: %10.3f\n", m_aspectRatio);
250 std::printf(" Loop closing threshold: %10.3f\n", m_loopThreshold);
251 if (m_sortContours) {
252 std::cout << " Sort contours.\n";
253 } else {
254 std::cout << " Do not sort contours.\n";
255 }
256 if (m_checkCrossings) {
257 std::cout << " Check for crossings.\n";
258 } else {
259 std::cout << " Do not check for crossings.\n";
260 }
261 if (markers) {
262 std::cout << " Mark isochron points.\n";
263 } else {
264 std::cout << " Draw isochron lines.\n";
265 }
266 }
267
268 // Loop over the equal time contours.
269 TGraph graph;
270 graph.SetLineColor(kGray + 2);
271 graph.SetMarkerColor(kGray + 2);
272 graph.SetLineStyle(m_lineStyle);
273 graph.SetMarkerStyle(m_markerStyle);
274 const double colRange = double(gStyle->GetNumberOfColors()) / nContours;
275 for (unsigned int ic = 0; ic < nContours; ++ic) {
276 if (colour) {
277 const auto col = gStyle->GetColorPalette(int((ic + 0.99) * colRange));
278 graph.SetLineColor(col);
279 graph.SetMarkerColor(col);
280 }
281 for (int stat : allStats) {
282 std::vector<std::pair<std::array<double, 4>, unsigned int> > contour;
283 // Loop over the drift lines, picking up the points when OK.
284 for (unsigned int k = 0; k < nDriftLines; ++k) {
285 const auto& dl = driftLines[k];
286 // Reject any undesirable combinations.
287 if (statusCodes[k] != stat || ic >= dl.size()) continue;
288 // Add the point to the contour line.
289 std::array<double, 4> point = {dl[ic][0], dl[ic][1], dl[ic][2], 0.};
290 contour.push_back(std::make_pair(point, k));
291 }
292 // Skip the plot of this contour if there are no points.
293 if (contour.empty()) continue;
294 bool circle = false;
295 // If requested, sort the points on the contour line.
296 if (m_sortContours && !markers) SortContour(contour, circle);
297 // Plot this contour.
298 if (markers) {
299 // Simply mark the contours if this was requested.
300 std::vector<double> xp;
301 std::vector<double> yp;
302 std::vector<double> zp;
303 for (const auto& point : contour) {
304 const double x = point.first[0];
305 const double y = point.first[1];
306 const double z = point.first[2];
307 xp.push_back(m_proj[0][0] * x + m_proj[1][0] * y + z * m_plane[0]);
308 yp.push_back(m_proj[0][1] * x + m_proj[1][1] * y + z * m_plane[1]);
309 zp.push_back(m_proj[0][2] * x + m_proj[1][2] * y + z * m_plane[2]);
310 }
311 graph.DrawGraph(xp.size(), xp.data(), yp.data(), "Psame");
312 continue;
313 }
314 // Regular plotting.
315 const double tolx = (m_xMaxPlot - m_xMinPlot) * m_connectionThreshold;
316 const double toly = (m_yMaxPlot - m_yMinPlot) * m_connectionThreshold;
317 // Flag to keep track if the segment is interrupted by a drift line
318 // or if it is too long.
319 bool gap = false;
320 // Coordinates to be plotted.
321 std::vector<double> xp;
322 std::vector<double> yp;
323 std::vector<double> zp;
324 const unsigned int nP = contour.size();
325 for (unsigned int i = 0; i < nP; ++i) {
326 gap = false;
327 const auto x0 = contour[i].first[0];
328 const auto y0 = contour[i].first[1];
329 const auto z0 = contour[i].first[2];
330 xp.push_back(m_proj[0][0] * x0 + m_proj[1][0] * y0 + z0 * m_plane[0]);
331 yp.push_back(m_proj[0][1] * x0 + m_proj[1][1] * y0 + z0 * m_plane[1]);
332 zp.push_back(m_proj[0][2] * x0 + m_proj[1][2] * y0 + z0 * m_plane[2]);
333 if (i == nP - 1) break;
334 const auto x1 = contour[i + 1].first[0];
335 const auto y1 = contour[i + 1].first[1];
336 // Reject contour segments which are long compared with AREA.
337 if (fabs(x1 - x0) > tolx || fabs(y1 - y0) > toly) gap = true;
338 // Get the indices of the drift lines corresponding
339 // to these two points on the contour line.
340 const auto i0 = contour[i].second;
341 const auto i1 = contour[i + 1].second;
342 // Set the BREAK flag if it crosses some stored drift line segment.
343 if (m_checkCrossings && !gap) {
344 for (unsigned int k = 0; k < nDriftLines; ++k) {
345 const auto& dl = driftLines[k];
346 for (unsigned int jc = 0; jc < dl.size(); ++jc) {
347 if ((i0 == k || i1 == k) && (jc == ic || jc + 1 == ic)) {
348 continue;
349 }
350 const auto& p0 = dl[jc];
351 const auto& p1 = jc == dl.size() - 1 ? endPoints[k] : dl[jc + 1];
352 if (Crossing(p0[0], p0[1], p1[0], p1[1], x0, y0, x1, y1)) {
353 gap = true;
354 break;
355 }
356 }
357 if (gap) break;
358 if ((i0 == k || i1 == k) && ic == 0) continue;
359 const auto& p0 = startPoints[k];
360 if (Crossing(p0[0], p0[1], dl[0][0], dl[0][1],
361 x0, y0, x1, y1)) {
362 gap = true;
363 break;
364 }
365 }
366 }
367 // If there has been a break, plot what we have already.
368 if (gap) {
369 if (xp.size() > 1) {
370 // Plot line.
371 graph.DrawGraph(xp.size(), xp.data(), yp.data(), "Lsame");
372 } else {
373 // Plot markers.
374 graph.DrawGraph(xp.size(), xp.data(), yp.data(), "Psame");
375 }
376 xp.clear();
377 yp.clear();
378 zp.clear();
379 }
380 }
381 // Plot the remainder.
382 if (xp.size() > 1) {
383 graph.DrawGraph(xp.size(), xp.data(), yp.data(), "Lsame");
384 } else if (!xp.empty()) {
385 graph.DrawGraph(xp.size(), xp.data(), yp.data(), "Psame");
386 }
387 }
388 }
389
390 gPad->Update();
391 if (!plotDriftLines) return;
392
393 graph.SetLineStyle(1);
394 if (m_particle == Particle::Electron) {
395 graph.SetLineColor(kOrange - 3);
396 } else {
397 graph.SetLineColor(kRed + 1);
398 }
399 for (unsigned int i = 0; i < nDriftLines; ++i) {
400 std::vector<double> xp;
401 std::vector<double> yp;
402 const double x0 = startPoints[i][0];
403 const double y0 = startPoints[i][1];
404 const double z0 = startPoints[i][2];
405 xp.push_back(m_proj[0][0] * x0 + m_proj[1][0] * y0 + z0 * m_plane[0]);
406 yp.push_back(m_proj[0][1] * x0 + m_proj[1][1] * y0 + z0 * m_plane[1]);
407 for (const auto& point : driftLines[i]) {
408 const double x = point[0];
409 const double y = point[1];
410 const double z = point[2];
411 xp.push_back(m_proj[0][0] * x + m_proj[1][0] * y + z * m_plane[0]);
412 yp.push_back(m_proj[0][1] * x + m_proj[1][1] * y + z * m_plane[1]);
413 }
414 const double x1 = endPoints[i][0];
415 const double y1 = endPoints[i][1];
416 const double z1 = endPoints[i][2];
417 xp.push_back(m_proj[0][0] * x1 + m_proj[1][0] * y1 + z1 * m_plane[0]);
418 yp.push_back(m_proj[0][1] * x1 + m_proj[1][1] * y1 + z1 * m_plane[1]);
419 graph.DrawGraph(xp.size(), xp.data(), yp.data(), "Lsame");
420 }
421}
422
423void ViewIsochrons::ComputeDriftLines(const double tstep,
424 const std::vector<std::array<double, 3> >& points,
425 std::vector<std::vector<std::array<double, 3> > >& driftLines,
426 std::vector<std::array<double, 3> >& startPoints,
427 std::vector<std::array<double, 3> >& endPoints,
428 std::vector<int>& statusCodes, const bool rev) {
429
430 DriftLineRKF drift;
431 Sensor sensor;
432 if (m_sensor) {
433 drift.SetSensor(m_sensor);
434 } else {
435 sensor.AddComponent(m_component);
436 drift.SetSensor(&sensor);
437 }
438 const double lx = 0.1 * fabs(m_xMaxPlot - m_xMinPlot);
439 const double ly = 0.1 * fabs(m_yMaxPlot - m_yMinPlot);
440 drift.SetMaximumStepSize(std::min(lx, ly));
441 drift.EnableSignalCalculation(false);
442 for (const auto& point : points) {
443 if (m_particle == Particle::Electron) {
444 if (m_positive) {
445 drift.DriftPositron(point[0], point[1], point[2], 0.);
446 } else {
447 drift.DriftElectron(point[0], point[1], point[2], 0.);
448 }
449 } else {
450 if (m_positive) {
451 drift.DriftIon(point[0], point[1], point[2], 0.);
452 } else {
453 drift.DriftNegativeIon(point[0], point[1], point[2], 0.);
454 }
455 }
456 const unsigned int nu = drift.GetNumberOfDriftLinePoints();
457 // Check that the drift line has enough points.
458 if (nu < 3) continue;
459 int status = 0;
460 double xf = 0., yf = 0., zf = 0., tf = 0.;
461 drift.GetEndPoint(xf, yf, zf, tf, status);
462 // Find the number of points to be stored.
463 const unsigned int nSteps = static_cast<unsigned int>(tf / tstep);
464 if (nSteps == 0) continue;
465 std::vector<double> xu(nu, 0.);
466 std::vector<double> yu(nu, 0.);
467 std::vector<double> zu(nu, 0.);
468 std::vector<double> tu(nu, 0.);
469 for (unsigned int i = 0; i < nu; ++i) {
470 drift.GetDriftLinePoint(i, xu[i], yu[i], zu[i], tu[i]);
471 }
472 if (rev) {
473 for (auto& t : tu) t = tf - t;
474 std::reverse(std::begin(xu), std::end(xu));
475 std::reverse(std::begin(yu), std::end(yu));
476 std::reverse(std::begin(zu), std::end(zu));
477 std::reverse(std::begin(tu), std::end(tu));
478 }
479 std::vector<std::array<double, 3> > tab;
480 // Interpolate at regular time intervals.
481 for (unsigned int i = 0; i < nSteps; ++i) {
482 const double t = (i + 1) * tstep;
483 // tab.push_back(PLACO3(Interpolate(xu, tu, t),
484 // Interpolate(yu, tu, t),
485 // Interpolate(zu, tu, t)));
486 std::array<double, 3> step = {Interpolate(xu, tu, t),
487 Interpolate(yu, tu, t),
488 Interpolate(zu, tu, t)};
489 tab.push_back(step);
490 }
491 driftLines.push_back(std::move(tab));
492 std::array<double, 3> start = {xu[0], yu[0], zu[0]};
493 std::array<double, 3> end = {xu[nu - 1], yu[nu - 1], zu[nu - 1]};
494 startPoints.push_back(std::move(start));
495 endPoints.push_back(std::move(end));
496 // Store the drift line return code.
497 if (rev) {
498 statusCodes.push_back(status);
499 } else {
500 statusCodes.push_back(0);
501 }
502 }
503}
504
505bool ViewIsochrons::SetPlotLimits() {
506
507 if (m_userPlotLimits) return true;
508 double xmin = 0., ymin = 0., xmax = 0., ymax = 0.;
509 if (m_userBox) {
510 if (PlotLimitsFromUserBox(xmin, ymin, xmax, ymax)) {
511 m_xMinPlot = xmin;
512 m_xMaxPlot = xmax;
513 m_yMinPlot = ymin;
514 m_yMaxPlot = ymax;
515 return true;
516 }
517 }
518 // Try to get the area/bounding box from the sensor/component.
519 bool ok = false;
520 if (m_sensor) {
521 ok = PlotLimits(m_sensor, xmin, ymin, xmax, ymax);
522 } else {
523 ok = PlotLimits(m_component, xmin, ymin, xmax, ymax);
524 }
525 if (ok) {
526 m_xMinPlot = xmin;
527 m_xMaxPlot = xmax;
528 m_yMinPlot = ymin;
529 m_yMaxPlot = ymax;
530 }
531 return ok;
532}
533
534void ViewIsochrons::SortContour(
535 std::vector<std::pair<std::array<double, 4>, unsigned int> >& contour,
536 bool& circle) {
537
538 if (contour.size() < 2) return;
539 // First compute the centre of gravity.
540 double xcog = 0.;
541 double ycog = 0.;
542 for (const auto& point : contour) {
543 xcog += point.first[0];
544 ycog += point.first[1];
545 }
546 const double scale = 1. / contour.size();
547 xcog *= scale;
548 ycog *= scale;
549 // Compute angles wrt to the centre of gravity and principal axes.
550 double sxx = 0.;
551 double sxy = 0.;
552 double syy = 0.;
553 for (const auto& point : contour) {
554 const double dx = point.first[0] - xcog;
555 const double dy = point.first[1] - ycog;
556 sxx += dx * dx;
557 sxy += dx * dy;
558 syy += dy * dy;
559 }
560 const double theta = 0.5 * atan2(2 * sxy, sxx - syy);
561 const double ct = cos(theta);
562 const double st = sin(theta);
563 // Evaluate dispersions around the principal axes.
564 sxx = 0.;
565 syy = 0.;
566 for (const auto& point : contour) {
567 const double dx = point.first[0] - xcog;
568 const double dy = point.first[1] - ycog;
569 sxx += fabs(+ct * dx + st * dy);
570 syy += fabs(-st * dx + ct * dy);
571 }
572 // Decide whether this is more linear or more circular.
573 if (fabs(sxx) > m_aspectRatio * fabs(syy) ||
574 fabs(syy) > m_aspectRatio * fabs(sxx)) {
575 circle = false;
576 } else {
577 circle = true;
578 }
579 // Set a sorting coordinate accordingly.
580 for (auto& point : contour) {
581 const double dx = point.first[0] - xcog;
582 const double dy = point.first[1] - ycog;
583 point.first[3] = circle ? atan2(dy, dx) : ct * dx + st * dy;
584 }
585 // Sort the points.
586 std::sort(contour.begin(), contour.end(),
587 [](const std::pair<std::array<double, 4>, int>& p1,
588 const std::pair<std::array<double, 4>, int>& p2) {
589 return p1.first[3] < p2.first[3]; }
590 );
591 if (!circle) return;
592 // For circles, perhaps add the first point to the end of the list.
593 // Compute breakpoint, total distance and maximum distance.
594 double dsum = 0.;
595 double dmax = -1.;
596 unsigned int imax = 0;
597 const unsigned int nPoints = contour.size();
598 for (unsigned int j = 0; j < nPoints; ++j) {
599 const auto& p1 = contour[j];
600 const auto& p0 = j > 0 ? contour[j - 1] : contour.back();
601 const double dx = p1.first[0] - p0.first[0];
602 const double dy = p1.first[1] - p0.first[1];
603 const double d = sqrt(dx * dx + dy * dy);
604 if (j > 0) dsum += d;
605 if (dmax < d) {
606 dmax = d;
607 imax = j;
608 }
609 }
610 if (dmax < m_loopThreshold * dsum) {
611 // If a true loop, close it.
612 contour.push_back(contour[0]);
613 } else {
614 circle = false;
615 if (imax > 0) {
616 // Shift the points to make a line.
617 std::rotate(contour.begin(), contour.begin() + imax, contour.end());
618 }
619 }
620}
621
622}
Abstract base class for components.
Definition: Component.hh:13
void EnableSignalCalculation(const bool on=true)
Switch on/off calculation of induced currents (default: enabled).
Definition: DriftLineRKF.hh:33
bool DriftNegativeIon(const double x0, const double y0, const double z0, const double t0)
bool DriftIon(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an ion with a given starting point.
void SetSensor(Sensor *s)
Set the sensor.
Definition: DriftLineRKF.cc:55
unsigned int GetNumberOfDriftLinePoints() const
Get the number of points of the most recent drift line.
Definition: DriftLineRKF.hh:90
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an electron with a given starting point.
bool DriftPositron(const double x0, const double y0, const double z0, const double t0)
void GetDriftLinePoint(const unsigned int i, double &x, double &y, double &z, double &t) const
Get the coordinates and time of a point along the most recent drift line.
void GetEndPoint(double &x, double &y, double &z, double &t, int &st) const
Get the end point and status flag of the most recent drift line.
void SetMaximumStepSize(const double ms)
Set the maximum step size that is allowed. By default, there is no limit.
Definition: DriftLineRKF.cc:72
void AddComponent(Component *comp)
Add a component.
Definition: Sensor.cc:316
Base class for visualization classes.
Definition: ViewBase.hh:18
std::array< double, 4 > m_plane
Definition: ViewBase.hh:99
std::string LabelY()
Definition: ViewBase.cc:483
std::string LabelX()
Definition: ViewBase.cc:420
std::string m_className
Definition: ViewBase.hh:78
std::array< std::array< double, 3 >, 3 > m_proj
Definition: ViewBase.hh:96
TPad * GetCanvas()
Retrieve the canvas.
Definition: ViewBase.cc:74
bool PlotLimits(Sensor *sensor, double &xmin, double &ymin, double &xmax, double &ymax) const
Definition: ViewBase.cc:608
bool PlotLimitsFromUserBox(double &xmin, double &ymin, double &xmax, double &ymax) const
Definition: ViewBase.cc:663
void SetComponent(Component *c)
Set the component.
void PlotIsochrons(const double tstep, const std::vector< std::array< double, 3 > > &points, const bool rev=false, const bool colour=false, const bool markers=false, const bool plotDriftLines=true)
void SetConnectionThreshold(const double thr)
void SetAspectRatioSwitch(const double ar)
ViewIsochrons()
Constructor.
void SetLoopThreshold(const double thr)
void SetSensor(Sensor *s)
Set the sensor.
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314