Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewMedium.cc
Go to the documentation of this file.
1#include <cmath>
2#include <iostream>
3#include <string>
4#include <algorithm>
5#include <limits>
6
7#include <TH1F.h>
8#include <TAxis.h>
9#include <TGraph.h>
10#include <TLatex.h>
11
13#include "Garfield/Medium.hh"
14#include "Garfield/Plotting.hh"
16
17namespace {
18
19int FindIndex(const std::vector<double>& fields, const double field,
20 const double eps) {
21
22 if (fields.empty()) return -1;
23 const int n = fields.size();
24 for (int i = 0; i < n; ++i) {
25 const double sum = fabs(fields[i]) + fabs(field);
26 const double tol = std::max(eps * sum, Garfield::Small);
27 if (fabs(fields[i] - field) < tol) return i;
28 }
29 return -1;
30}
31
32bool NonZero(const std::vector<double>& v) {
33
34 constexpr double tol = 1.e-10;
35 return std::any_of(v.cbegin(), v.cend(), [](double x){ return fabs(x) > tol; });
36}
37
38}
39
40namespace Garfield {
41
42ViewMedium::ViewMedium() : ViewBase("ViewMedium") {}
43
45 if (!m) {
46 std::cerr << m_className << "::SetMedium: Null pointer.\n";
47 return;
48 }
49
50 m_medium = m;
51}
52
53void ViewMedium::SetRangeE(const double emin, const double emax,
54 const bool logscale) {
55 if (emin >= emax || emin < 0.) {
56 std::cerr << m_className << "::SetRangeE: Incorrect range.\n";
57 return;
58 }
59
60 m_eMin = emin;
61 m_eMax = emax;
62 m_logE = logscale;
63}
64
65void ViewMedium::SetRangeB(const double bmin, const double bmax, const bool logscale) {
66 if (bmin >= bmax || bmin < 0.) {
67 std::cerr << m_className << "::SetRangeB: Incorrect range.\n";
68 return;
69 }
70
71 m_bMin = bmin;
72 m_bMax = bmax;
73 m_logB = logscale;
74}
75
76void ViewMedium::SetRangeA(const double amin, const double amax, const bool logscale) {
77 if (amin >= amax || amin < 0.) {
78 std::cerr << m_className << "::SetRangeA: Incorrect range.\n";
79 return;
80 }
81
82 m_aMin = amin;
83 m_aMax = amax;
84 m_logA = logscale;
85}
86
87void ViewMedium::SetRangeY(const double ymin, const double ymax,
88 const bool logscale) {
89 if (ymin >= ymax || ymin < 0.) {
90 std::cerr << m_className << "::SetRangeY: Incorrect range.\n";
91 return;
92 }
93
94 m_yMin = ymin;
95 m_yMax = ymax;
96 m_logY = logscale;
97}
98
100 std::cerr << m_className << "::PlotElectronCrossSections: Not implemented.\n";
101}
102
103void ViewMedium::Draw() {
104
105 if (m_yPlot.empty()) return;
106 auto canvas = GetCanvas();
107 canvas->cd();
108
109 // Find the x and y ranges.
110 const double xmin = m_xPlot.front();
111 const double xmax = m_xPlot.back();
112 double ymin = m_yMin;
113 double ymax = m_yMax;
114 if (m_autoRangeY) {
115 ymin = std::numeric_limits<double>::max();
116 ymax = -std::numeric_limits<double>::max();
117 for (const auto& plot : m_yPlot) {
118 ymin = std::min(ymin, *std::min_element(plot.cbegin(), plot.cend()));
119 ymax = std::max(ymax, *std::max_element(plot.cbegin(), plot.cend()));
120 }
121 const double dy = 0.05 * fabs(ymax - ymin);
122 ymin = std::max(0., ymin - dy);
123 ymax += dy;
124 }
125 auto frame = canvas->DrawFrame(xmin, ymin, xmax, ymax);
126 // Set the axis labels.
127 if (m_xaxis == Axis::E) {
128 frame->GetXaxis()->SetTitle("electric field [V/cm]");
129 } else if (m_xaxis == Axis::B) {
130 frame->GetXaxis()->SetTitle("magnetic field [T]");
131 } else if (m_xaxis == Axis::Angle) {
132 frame->GetXaxis()->SetTitle("angle between #bf{E} and #bf{B} [rad]");
133 } else {
134 std::cerr << m_className << "::Draw: Invalid x axis.\n";
135 return;
136 }
137 auto yaxis = frame->GetYaxis();
138 switch (m_par[0]) {
139 case Parameter::VelocityE:
140 case Parameter::VelocityB:
141 case Parameter::VelocityExB:
142 yaxis->SetTitle("drift velocity [cm/ns]");
143 canvas->SetTitle("Drift velocity");
144 break;
145 case Parameter::LongitudinalDiffusion:
146 case Parameter::TransverseDiffusion:
147 yaxis->SetTitle("diffusion coefficient [#kern[-0.1]{#sqrt{cm}}#kern[0.1]{]}");
148 canvas->SetTitle("Diffusion");
149 break;
150 case Parameter::Townsend:
151 case Parameter::Attachment:
152 if (std::all_of(m_par.cbegin(), m_par.cend(), [](Parameter p) {
153 return p == Parameter::Townsend; })) {
154 yaxis->SetTitle("#it{#alpha} [1/cm]");
155 canvas->SetTitle("Multiplication");
156 } else if (std::all_of(m_par.cbegin(), m_par.cend(), [](Parameter p) {
157 return p == Parameter::Attachment; })) {
158 yaxis->SetTitle("#it{#eta} [1/cm]");
159 canvas->SetTitle("Attachment");
160 } else {
161 yaxis->SetTitle("#it{#alpha}, #it{#eta} [1/cm]");
162 canvas->SetTitle("Multiplication and attachment");
163 }
164 break;
165 case Parameter::LorentzAngle:
166 yaxis->SetTitle( "Angle between #bf{v} and #bf{E} [rad]");
167 canvas->SetTitle("Lorentz angle");
168 break;
169 default:
170 break;
171 }
172 yaxis->SetTitleOffset(0);
173 gPad->SetLeftMargin(0.15);
174 canvas->Update();
175
176 const unsigned int nPlots = m_yPlot.size();
177 // Set colours.
178 std::vector<short> cols(nPlots, 0);
179 if (m_colours.empty()) {
180 auto it = std::find(m_q.cbegin(), m_q.cend(), Charge::Electron);
181 if (it != m_q.cend()) {
182 cols[std::distance(m_q.cbegin(), it)] = kOrange - 3;
183 }
184 it = std::find(std::next(it), m_q.cend(), Charge::Electron);
185 if (it != m_q.cend()) {
186 cols[std::distance(m_q.cbegin(), it)] = kGreen + 3;
187 }
188 it = std::find(m_q.cbegin(), m_q.cend(), Charge::Hole);
189 if (it != m_q.cend()) {
190 cols[std::distance(m_q.cbegin(), it)] = kRed + 1;
191 }
192 it = std::find(m_q.cbegin(), m_q.cend(), Charge::Ion);
193 if (it != m_q.cend()) {
194 cols[std::distance(m_q.cbegin(), it)] = kRed + 1;
195 }
196 } else {
197 for (unsigned int i = 0; i < nPlots; ++i) {
198 cols[i] = m_colours[i % m_colours.size()];
199 }
200 }
201 // Set legend.
202 std::vector<std::string> labels = m_labels;
203 labels.resize(nPlots, "");
204 if (nPlots > 1 && m_labels.empty()) {
205 bool allEqual = std::equal(m_q.begin() + 1, m_q.end(), m_q.begin());
206 if (!allEqual) {
207 for (unsigned int i = 0; i < nPlots; ++i) {
208 if (m_q[i] == Charge::Electron) {
209 labels[i] = "electrons";
210 } else if (m_q[i] == Charge::Hole) {
211 labels[i] = "holes";
212 } else {
213 labels[i] = "ions";
214 }
215 }
216 }
217 allEqual = std::equal(m_par.begin() + 1, m_par.end(), m_par.begin());
218 if (!allEqual) {
219 for (unsigned int i = 0; i < nPlots; ++i) {
220 if (!labels[i].empty()) labels[i] += ", ";
221 switch (m_par[i]) {
222 case Parameter::VelocityE:
223 labels[i] += "#it{v}_{#it{E}}";
224 break;
225 case Parameter::VelocityB:
226 labels[i] += "#it{v}_{#it{B}#kern[0.1]{t}}";
227 break;
228 case Parameter::VelocityExB:
229 labels[i] += "#it{v}_{#it{E}#kern[0.05]{#times}#it{B}}";
230 break;
231 case Parameter::LongitudinalDiffusion:
232 labels[i] += "#it{D}_{L}";
233 break;
234 case Parameter::TransverseDiffusion:
235 labels[i] += "#it{D}_{T}";
236 break;
237 case Parameter::Townsend:
238 labels[i] += "#alpha";
239 break;
240 case Parameter::Attachment:
241 labels[i] += "#eta";
242 break;
243 case Parameter::LorentzAngle:
244 labels[i] += "Lorentz angle";
245 }
246 }
247 }
248 }
249 TGraph graph;
250 TLatex latex;
251 double xsize = 0., ysize = 0.;
252 for (const auto& label : labels) {
253 latex.SetText(0, 0, label.c_str());
254 xsize = std::max(xsize, latex.GetXsize());
255 ysize = std::max(ysize, latex.GetYsize());
256 }
257 // Convert to NDC.
258 xsize /= (gPad->GetX2() - gPad->GetX1());
259 ysize /= (gPad->GetY2() - gPad->GetY1());
260
261 const double lm = gPad->GetLeftMargin();
262 const double rm = 1. - gPad->GetRightMargin();
263 const double tm = 1. - gPad->GetTopMargin();
264 const double bm = gPad->GetBottomMargin();
265 double xLabel = lm + 0.1 * (rm - lm);
266 if (m_par[0] == Parameter::LongitudinalDiffusion ||
267 m_par[0] == Parameter::TransverseDiffusion) {
268 xLabel = rm - 0.1 * (rm - lm) - xsize;
269 }
270 double yLabel = tm - 0.1 * (tm - bm);
271 const double colrange = gStyle->GetNumberOfColors() / double(nPlots);
272 for (unsigned int i = 0; i < nPlots; ++i) {
273 int col = cols[i] > 0 ? cols[i] : gStyle->GetColorPalette(i * colrange);
274 graph.SetLineColor(col);
275 graph.SetMarkerColor(col);
276 graph.DrawGraph(m_xPlot.size(),
277 m_xPlot.data(), m_yPlot[i].data(), "L");
278 graph.SetMarkerStyle(20 + i);
279 if (!m_xGraph[i].empty()) {
280 graph.DrawGraph(m_xGraph[i].size(),
281 m_xGraph[i].data(), m_yGraph[i].data(), "P");
282
283 }
284 if (labels[i].empty()) continue;
285 latex.SetTextColor(col);
286 latex.DrawLatexNDC(xLabel, yLabel, labels[i].c_str());
287 yLabel -= 1.5 * ysize;
288 }
289 if (m_logX && xmin > 0.) {
290 canvas->SetLogx(1);
291 } else {
292 canvas->SetLogx(0);
293 }
294 if (m_logY && ymin > 0. && !m_autoRangeY) {
295 canvas->SetLogy(1);
296 } else {
297 canvas->SetLogy(0);
298 }
299 gPad->Update();
300}
301
302void ViewMedium::ResetY() {
303 m_yPlot.clear();
304 m_par.clear();
305 m_q.clear();
306 m_xGraph.clear();
307 m_yGraph.clear();
308}
309
310void ViewMedium::ResetX(const Axis xaxis) {
311
312 double xmin = 0., xmax = 0.;
313 bool logx = false;
314 if (xaxis == Axis::E) {
315 xmin = m_eMin;
316 xmax = m_eMax;
317 logx = m_logE;
318 } else if (xaxis == Axis::B) {
319 xmin = m_bMin;
320 xmax = m_bMax;
321 logx = m_logB;
322 } else if (xaxis == Axis::Angle) {
323 xmin = m_aMin;
324 xmax = m_aMax;
325 logx = m_logA;
326 } else {
327 m_xaxis = Axis::None;
328 return;
329 }
330
331 if (m_autoRangeX) {
332 // Get the field grid.
333 std::vector<double> efields;
334 std::vector<double> bfields;
335 std::vector<double> bangles;
336 m_medium->GetFieldGrid(efields, bfields, bangles);
337 if (xaxis == Axis::E && !efields.empty()) {
338 if (efields.front() > 0. && efields.back() > 10. * efields.front()) {
339 logx = true;
340 xmin = efields.front() / 1.5;
341 xmax = efields.back() * 1.5;
342 } else {
343 logx = false;
344 const double dx = 0.05 * fabs(efields.back() - efields.front());
345 xmin = std::max(0., efields.front() - dx);
346 xmax = efields.back() + dx;
347 }
348 } else if (xaxis == Axis::B && !bfields.empty()) {
349 logx = false;
350 const double dx = 0.05 * fabs(bfields.back() - bfields.front());
351 xmin = std::max(0., bfields.front() - dx);
352 xmax = bfields.back() + dx;
353 } else if (xaxis == Axis::Angle && !bangles.empty()) {
354 logx = false;
355 const double dx = 0.05 * fabs(bangles.back() - bangles.front());
356 xmin = std::max(0., bangles.front() - dx);
357 xmax = std::min(TwoPi, bangles.back() + dx);
358 }
359 }
360
361 constexpr unsigned int nX = 1000;
362 m_xPlot.assign(nX, 0.);
363 if (logx) {
364 m_xPlot[0] = xmin;
365 const double r = pow(xmax / xmin, 1. / (nX - 1));
366 for (unsigned int i = 1; i < nX; ++i) {
367 m_xPlot[i] = m_xPlot[i - 1] * r;
368 }
369 } else {
370 const double dx = (xmax - xmin) / (nX - 1);
371 for (unsigned int i = 0; i < nX; ++i) {
372 m_xPlot[i] = xmin + i * dx;
373 }
374 }
375 m_logX = logx;
376 m_xaxis = xaxis;
377}
378
379void ViewMedium::PlotDiffusion(const Axis xaxis, const Charge charge,
380 const bool same) {
381
382 if (!m_medium) {
383 std::cerr << m_className << "::PlotDiffusion: Medium is not defined.\n";
384 return;
385 }
386 if (xaxis != m_xaxis) {
387 ResetX(xaxis);
388 ResetY();
389 } else if (!same) {
390 ResetY();
391 } else if (!m_par.empty()) {
392 if (m_par[0] != Parameter::TransverseDiffusion &&
393 m_par[0] != Parameter::LongitudinalDiffusion) {
394 ResetY();
395 }
396 }
397 const unsigned int nX = m_xPlot.size();
398 std::array<std::vector<double>, 2> ypl;
399 for (unsigned int i = 0; i < 2; ++i) ypl[i].assign(nX, 0.);
400
401 double ex = m_efield;
402 double ctheta = cos(m_angle);
403 double stheta = sin(m_angle);
404 double bx = m_bfield * ctheta;
405 double by = m_bfield * stheta;
406 for (unsigned int i = 0; i < nX; ++i) {
407 if (xaxis == Axis::E) {
408 ex = m_xPlot[i];
409 } else if (xaxis == Axis::B) {
410 bx = m_xPlot[i] * ctheta;
411 by = m_xPlot[i] * stheta;
412 } else {
413 bx = m_bfield * cos(m_xPlot[i]);
414 by = m_bfield * sin(m_xPlot[i]);
415 }
416 double dl = 0., dt = 0.;
417 if (charge == Charge::Electron) {
418 if (!m_medium->ElectronDiffusion(ex, 0, 0, bx, by, 0, dl, dt)) continue;
419 } else if (charge == Charge::Hole) {
420 if (!m_medium->HoleDiffusion(ex, 0, 0, bx, by, 0, dl, dt)) continue;
421 } else {
422 if (!m_medium->IonDiffusion(ex, 0, 0, bx, by, 0, dl, dt)) continue;
423 }
424 ypl[0][i] = dl;
425 ypl[1][i] = dt;
426 }
427
428 std::array<std::vector<double>, 2> xgr;
429 std::array<std::vector<double>, 2> ygr;
430
431 std::array<std::vector<double>, 3> grid;
432 int ie = 0, ib = 0, ia = 0;
433 if (GetGrid(grid, ie, ib, ia, xaxis)) {
434 const auto nPoints = xaxis == Axis::E ? grid[0].size() :
435 xaxis == Axis::B ? grid[1].size() : grid[2].size();
436 for (unsigned int j = 0; j < nPoints; ++j) {
437 double x = 0., y = 0.;
438 if (xaxis == Axis::E) {
439 ie = j;
440 x = m_medium->UnScaleElectricField(grid[0][j]);
441 } else if (xaxis == Axis::B) {
442 ib = j;
443 x = grid[1][j];
444 } else if (xaxis == Axis::Angle) {
445 ia = j;
446 x = grid[2][j];
447 }
448 if (charge == Charge::Electron) {
449 if (m_medium->GetElectronTransverseDiffusion(ie, ib, ia, y)) {
450 xgr[1].push_back(x);
451 ygr[1].push_back(m_medium->ScaleDiffusion(y));
452 }
453 if (m_medium->GetElectronLongitudinalDiffusion(ie, ib, ia, y)) {
454 xgr[0].push_back(x);
455 ygr[0].push_back(m_medium->ScaleDiffusion(y));
456 }
457 } else if (charge == Charge::Hole) {
458 if (m_medium->GetHoleTransverseDiffusion(ie, ib, ia, y)) {
459 xgr[1].push_back(x);
460 ygr[1].push_back(m_medium->ScaleDiffusion(y));
461 }
462 if (m_medium->GetHoleLongitudinalDiffusion(ie, ib, ia, y)) {
463 xgr[0].push_back(x);
464 ygr[0].push_back(m_medium->ScaleDiffusion(y));
465 }
466 } else {
467 if (m_medium->GetIonTransverseDiffusion(ie, ib, ia, y)) {
468 xgr[1].push_back(x);
469 ygr[1].push_back(m_medium->ScaleDiffusion(y));
470 }
471 if (m_medium->GetIonLongitudinalDiffusion(ie, ib, ia, y)) {
472 xgr[0].push_back(x);
473 ygr[0].push_back(m_medium->ScaleDiffusion(y));
474 }
475 }
476 }
477 }
478
479 const std::array<Parameter, 2> pars = {Parameter::LongitudinalDiffusion,
480 Parameter::TransverseDiffusion};
481 for (unsigned int i = 0; i < 2; ++i) {
482 if (!NonZero(ypl[i])) continue;
483 m_yPlot.push_back(std::move(ypl[i]));
484 m_par.push_back(pars[i]);
485 m_q.push_back(charge);
486 m_xGraph.push_back(std::move(xgr[i]));
487 m_yGraph.push_back(std::move(ygr[i]));
488 }
489 Draw();
490}
491
492void ViewMedium::PlotVelocity(const Axis xaxis, const Charge charge,
493 const bool same) {
494
495 if (!m_medium) {
496 std::cerr << m_className << "::PlotVelocity: Medium is not defined.\n";
497 return;
498 }
499 if (xaxis != m_xaxis) {
500 ResetX(xaxis);
501 ResetY();
502 } else if (!same) {
503 ResetY();
504 } else if (!m_par.empty()) {
505 if (m_par[0] != Parameter::VelocityE &&
506 m_par[0] != Parameter::VelocityB &&
507 m_par[0] != Parameter::VelocityExB) {
508 ResetY();
509 }
510 }
511 const unsigned int nX = m_xPlot.size();
512 std::array<std::vector<double>, 3> ypl;
513 for (unsigned int i = 0; i < 3; ++i) ypl[i].assign(nX, 0.);
514
515 double e0 = m_efield;
516 double b0 = m_bfield;
517 double ctheta = cos(m_angle);
518 double stheta = sin(m_angle);
519
520 for (unsigned int i = 0; i < nX; ++i) {
521 if (xaxis == Axis::E) {
522 e0 = m_xPlot[i];
523 } else if (xaxis == Axis::B) {
524 b0 = m_xPlot[i];
525 } else {
526 ctheta = cos(m_xPlot[i]);
527 stheta = sin(m_xPlot[i]);
528 }
529 double vx = 0., vy = 0., vz = 0.;
530 if (charge == Charge::Electron) {
531 if (m_medium->ElectronVelocity(e0, 0, 0, b0 * ctheta, b0 * stheta, 0,
532 vx, vy, vz)) {
533 ypl[0][i] = fabs(vx);
534 ypl[1][i] = fabs(vy);
535 ypl[2][i] = fabs(vz);
536 }
537 } else if (charge == Charge::Hole) {
538 if (m_medium->HoleVelocity(e0, 0, 0, b0 * ctheta, b0 * stheta, 0,
539 vx, vy, vz)) {
540 ypl[0][i] = fabs(vx);
541 ypl[1][i] = fabs(vy);
542 ypl[2][i] = fabs(vz);
543 }
544 } else {
545 if (m_medium->IonVelocity(e0, 0, 0, b0 * ctheta, b0 * stheta, 0,
546 vx, vy, vz)) {
547 ypl[0][i] = fabs(vx);
548 }
549 }
550 }
551 std::array<std::vector<double>, 3> xgr;
552 std::array<std::vector<double>, 3> ygr;
553
554 std::array<std::vector<double>, 3> grid;
555 int ie = 0, ib = 0, ia = 0;
556 if (GetGrid(grid, ie, ib, ia, xaxis)) {
557 const auto nPoints = xaxis == Axis::E ? grid[0].size() :
558 xaxis == Axis::B ? grid[1].size() : grid[2].size();
559 for (unsigned int j = 0; j < nPoints; ++j) {
560 double x = 0., y = 0.;
561 if (xaxis == Axis::E) {
562 ie = j;
563 x = m_medium->UnScaleElectricField(grid[0][j]);
564 } else if (xaxis == Axis::B) {
565 ib = j;
566 x = grid[1][j];
567 } else if (xaxis == Axis::Angle) {
568 ia = j;
569 x = grid[2][j];
570 }
571 if (charge == Charge::Electron) {
572 if (m_medium->GetElectronVelocityE(ie, ib, ia, y)) {
573 xgr[0].push_back(x);
574 ygr[0].push_back(m_medium->ScaleVelocity(y));
575 }
576 if (m_medium->GetElectronVelocityB(ie, ib, ia, y)) {
577 xgr[1].push_back(x);
578 ygr[1].push_back(fabs(m_medium->ScaleVelocity(y)));
579 }
580 if (m_medium->GetElectronVelocityExB(ie, ib, ia, y)) {
581 xgr[2].push_back(x);
582 ygr[2].push_back(fabs(m_medium->ScaleVelocity(y)));
583 }
584 } else if (charge == Charge::Hole) {
585 if (m_medium->GetHoleVelocityE(ie, ib, ia, y)) {
586 xgr[0].push_back(x);
587 ygr[0].push_back(m_medium->ScaleVelocity(y));
588 }
589 if (m_medium->GetHoleVelocityB(ie, ib, ia, y)) {
590 xgr[1].push_back(x);
591 ygr[1].push_back(fabs(m_medium->ScaleVelocity(y)));
592 }
593 if (m_medium->GetHoleVelocityExB(ie, ib, ia, y)) {
594 xgr[2].push_back(x);
595 ygr[2].push_back(fabs(m_medium->ScaleVelocity(y)));
596 }
597 } else {
598 if (m_medium->GetIonMobility(ie, ib, ia, y)) {
599 xgr[0].push_back(x);
600 ygr[0].push_back(y * m_medium->UnScaleElectricField(grid[0][ie]));
601 }
602 }
603 }
604 }
605
606 const std::array<Parameter, 3> pars = {Parameter::VelocityE,
607 Parameter::VelocityB,
608 Parameter::VelocityExB};
609 for (unsigned int i = 0; i < 3; ++i) {
610 if (!NonZero(ypl[i])) continue;
611 m_yPlot.push_back(std::move(ypl[i]));
612 m_par.push_back(pars[i]);
613 m_q.push_back(charge);
614 m_xGraph.push_back(std::move(xgr[i]));
615 m_yGraph.push_back(std::move(ygr[i]));
616 }
617 Draw();
618}
619
620void ViewMedium::Plot(const Axis xaxis, const Charge charge,
621 const Parameter par, const bool same) {
622
623 // Make sure the medium is set.
624 if (!m_medium) {
625 std::cerr << m_className << "::Plot: Medium is not defined.\n";
626 return;
627 }
628 if (xaxis != m_xaxis) {
629 ResetX(xaxis);
630 ResetY();
631 } else if (!same) {
632 ResetY();
633 } else if (!m_par.empty()) {
634 if (m_par[0] != Parameter::Townsend &&
635 m_par[0] != Parameter::Attachment) {
636 ResetY();
637 }
638 }
639
640 const unsigned int nX = m_xPlot.size();
641 std::vector<double> ypl(nX, 0.);
642 double ex = m_efield;
643 double ctheta = cos(m_angle);
644 double stheta = sin(m_angle);
645 double bx = m_bfield * ctheta;
646 double by = m_bfield * stheta;
647 for (unsigned int i = 0; i < nX; ++i) {
648 if (xaxis == Axis::E) {
649 ex = m_xPlot[i];
650 } else if (xaxis == Axis::B) {
651 bx = m_xPlot[i] * ctheta;
652 by = m_xPlot[i] * stheta;
653 } else {
654 bx = m_bfield * cos(m_xPlot[i]);
655 by = m_bfield * sin(m_xPlot[i]);
656 }
657 double y = 0.;
658 if (charge == Charge::Electron) {
659 if (par == Parameter::Townsend) {
660 if (!m_medium->ElectronTownsend(ex, 0, 0, bx, by, 0, y)) continue;
661 } else {
662 if (!m_medium->ElectronAttachment(ex, 0, 0, bx, by, 0, y)) continue;
663 y = std::abs(y);
664 }
665 } else {
666 if (par == Parameter::Townsend) {
667 if (!m_medium->HoleTownsend(ex, 0, 0, bx, by, 0, y)) continue;
668 } else {
669 if (!m_medium->HoleAttachment(ex, 0, 0, bx, by, 0, y)) continue;
670 y = std::abs(y);
671 }
672 }
673 ypl[i] = y;
674 }
675
676 std::vector<double> xgr;
677 std::vector<double> ygr;
678
679 std::array<std::vector<double>, 3> grid;
680 int ie = 0, ib = 0, ia = 0;
681 if (GetGrid(grid, ie, ib, ia, xaxis)) {
682 const auto nPoints = xaxis == Axis::E ? grid[0].size() :
683 xaxis == Axis::B ? grid[1].size() : grid[2].size();
684 for (unsigned int j = 0; j < nPoints; ++j) {
685 double x = 0., y = 0.;
686 if (xaxis == Axis::E) {
687 ie = j;
688 x = m_medium->UnScaleElectricField(grid[0][j]);
689 } else if (xaxis == Axis::B) {
690 ib = j;
691 x = grid[1][j];
692 } else if (xaxis == Axis::Angle) {
693 ia = j;
694 x = grid[2][j];
695 }
696 if (charge == Charge::Electron) {
697 if (par == Parameter::Townsend) {
698 if (m_medium->GetElectronTownsend(ie, ib, ia, y)) {
699 xgr.push_back(x);
700 ygr.push_back(m_medium->ScaleTownsend(exp(y)));
701 }
702 } else {
703 if (m_medium->GetElectronAttachment(ie, ib, ia, y)) {
704 xgr.push_back(x);
705 ygr.push_back(m_medium->ScaleAttachment(exp(y)));
706 }
707 }
708 } else if (charge == Charge::Hole) {
709 if (par == Parameter::Townsend) {
710 if (m_medium->GetHoleTownsend(ie, ib, ia, y)) {
711 xgr.push_back(x);
712 ygr.push_back(m_medium->ScaleTownsend(exp(y)));
713 }
714 } else {
715 if (m_medium->GetHoleAttachment(ie, ib, ia, y)) {
716 xgr.push_back(x);
717 ygr.push_back(m_medium->ScaleAttachment(exp(y)));
718 }
719 }
720 }
721 }
722 }
723
724 if (!NonZero(ypl)) return;
725 m_yPlot.push_back(std::move(ypl));
726 m_par.push_back(par);
727 m_q.push_back(charge);
728 m_xGraph.push_back(std::move(xgr));
729 m_yGraph.push_back(std::move(ygr));
730 Draw();
731}
732
733void ViewMedium::PlotLorentzAngle(const Axis xaxis, const Charge charge,
734 const bool same) {
735
736 // Make sure the medium is set.
737 if (!m_medium) {
738 std::cerr << m_className << "::PlotLorentzAngle: Medium is not defined.\n";
739 return;
740 }
741 if (xaxis != m_xaxis) {
742 ResetX(xaxis);
743 ResetY();
744 } else if (!same) {
745 ResetY();
746 } else if (!m_par.empty()) {
747 if (m_par[0] != Parameter::LorentzAngle) ResetY();
748 }
749
750 const unsigned int nX = m_xPlot.size();
751 std::vector<double> ypl(nX, 0.);
752 double ex = m_efield;
753 double ctheta = cos(m_angle);
754 double stheta = sin(m_angle);
755 double bx = m_bfield * ctheta;
756 double by = m_bfield * stheta;
757 for (unsigned int i = 0; i < nX; ++i) {
758 if (xaxis == Axis::E) {
759 ex = m_xPlot[i];
760 } else if (xaxis == Axis::B) {
761 bx = m_xPlot[i] * ctheta;
762 by = m_xPlot[i] * stheta;
763 } else {
764 bx = m_bfield * cos(m_xPlot[i]);
765 by = m_bfield * sin(m_xPlot[i]);
766 }
767 double y = 0.;
768 if (!m_medium->ElectronLorentzAngle(ex, 0, 0, bx, by, 0, y)) continue;
769 ypl[i] = y;
770 }
771
772 std::vector<double> xgr;
773 std::vector<double> ygr;
774
775 std::array<std::vector<double>, 3> grid;
776 int ie = 0, ib = 0, ia = 0;
777 if (GetGrid(grid, ie, ib, ia, xaxis)) {
778 const auto nPoints = xaxis == Axis::E ? grid[0].size() :
779 xaxis == Axis::B ? grid[1].size() : grid[2].size();
780 for (unsigned int j = 0; j < nPoints; ++j) {
781 double x = 0., y = 0.;
782 if (xaxis == Axis::E) {
783 ie = j;
784 x = m_medium->UnScaleElectricField(grid[0][j]);
785 } else if (xaxis == Axis::B) {
786 ib = j;
787 x = grid[1][j];
788 } else if (xaxis == Axis::Angle) {
789 ia = j;
790 x = grid[2][j];
791 }
792 if (m_medium->GetElectronLorentzAngle(ie, ib, ia, y)) {
793 xgr.push_back(x);
794 ygr.push_back(y);
795 }
796 }
797 }
798
799 if (!NonZero(ypl)) return;
800 m_yPlot.push_back(std::move(ypl));
801 m_par.push_back(Parameter::LorentzAngle);
802 m_q.push_back(charge);
803 m_xGraph.push_back(std::move(xgr));
804 m_yGraph.push_back(std::move(ygr));
805 Draw();
806}
807
808ViewMedium::Axis ViewMedium::GetAxis(const char xaxis) const {
809
810 if (std::toupper(xaxis) == 'E') {
811 return Axis::E;
812 } else if (std::toupper(xaxis) == 'B') {
813 return Axis::B;
814 } else if (std::toupper(xaxis) == 'A') {
815 return Axis::Angle;
816 }
817 return Axis::None;
818}
819
820bool ViewMedium::GetGrid(std::array<std::vector<double>, 3>& grid,
821 int& ie, int& ib, int& ia,
822 const Axis xaxis) const {
823
824 if (!m_medium) return false;
825 m_medium->GetFieldGrid(grid[0], grid[1], grid[2]);
826 if (grid[0].empty() || grid[1].empty() || grid[2].empty()) return false;
827 constexpr double eps = 1.e-3;
828 ie = FindIndex(grid[0], m_efield, eps);
829 ib = FindIndex(grid[1], m_bfield, eps);
830 ia = FindIndex(grid[2], m_angle, eps);
831 if (xaxis == Axis::E) {
832 if (ib < 0 || ia < 0) return false;
833 } else if (xaxis == Axis::B) {
834 if (ie < 0 || ia < 0) return false;
835 } else if (xaxis == Axis::Angle) {
836 if (ie < 0 || ib < 0) return false;
837 } else {
838 return false;
839 }
840 return true;
841}
842
843}
Abstract base class for media.
Definition: Medium.hh:13
virtual double UnScaleElectricField(const double e) const
Definition: Medium.hh:464
bool GetHoleLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dl)
Get an entry in the table of longitudinal diffusion coefficients.
Definition: Medium.hh:321
virtual double ScaleVelocity(const double v) const
Definition: Medium.hh:465
bool GetHoleTownsend(const size_t ie, const size_t ib, const size_t ia, double &alpha)
Get an entry in the table of Townsend coefficients.
Definition: Medium.hh:341
virtual bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:534
bool GetElectronLorentzAngle(const size_t ie, const size_t ib, const size_t ia, double &lor)
Get an entry in the table of Lorentz angles.
Definition: Medium.hh:279
virtual double ScaleAttachment(const double eta) const
Definition: Medium.hh:469
virtual double ScaleDiffusion(const double d) const
Definition: Medium.hh:466
bool GetHoleAttachment(const size_t ie, const size_t ib, const size_t ia, double &eta)
Get an entry in the table of attachment coefficients.
Definition: Medium.hh:351
virtual bool HoleVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:513
virtual bool ElectronLorentzAngle(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
Lorentz angle.
Definition: Medium.cc:429
bool GetIonTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dt)
Get an entry in the table of transverse diffusion coefficients.
Definition: Medium.hh:387
bool GetElectronTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dt)
Get an entry in the table of transverse diffusion coefficients.
Definition: Medium.hh:248
virtual bool ElectronVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:379
virtual bool IonVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:565
bool GetHoleTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dt)
Get an entry in the table of transverse diffusion coefficients.
Definition: Medium.hh:331
bool GetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along ExB.
Definition: Medium.hh:300
virtual double ScaleTownsend(const double alpha) const
Definition: Medium.hh:468
bool GetElectronAttachment(const size_t ie, const size_t ib, const size_t ia, double &eta)
Get an entry in the table of attachment coefficients.
Definition: Medium.hh:268
virtual bool ElectronDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Longitudinal and transverse diffusion coefficients [cm1/2].
Definition: Medium.cc:387
bool GetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along E.
Definition: Medium.hh:208
void GetFieldGrid(std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles)
Get the fields and E-B angles used in the transport tables.
Definition: Medium.cc:860
virtual bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:403
bool GetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along Btrans.
Definition: Medium.hh:310
bool GetElectronTownsend(const size_t ie, const size_t ib, const size_t ia, double &alpha)
Get an entry in the table of Townsend coefficients.
Definition: Medium.hh:258
bool GetIonMobility(const size_t ie, const size_t ib, const size_t ia, double &mu)
Get an entry in the table of ion mobilities.
Definition: Medium.hh:366
virtual bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:416
bool GetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along ExB.
Definition: Medium.hh:218
bool GetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along Btrans.
Definition: Medium.hh:228
bool GetIonLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dl)
Get an entry in the table of longitudinal diffusion coefficients.
Definition: Medium.hh:377
virtual bool HoleDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Longitudinal and transverse diffusion coefficients [cm1/2].
Definition: Medium.cc:521
bool GetElectronLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dl)
Get an entry in the table of longitudinal diffusion coefficients.
Definition: Medium.hh:238
bool GetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along E.
Definition: Medium.hh:290
virtual bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:547
virtual bool IonDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Longitudinal and transverse diffusion coefficients [cm1/2].
Definition: Medium.cc:600
Base class for visualization classes.
Definition: ViewBase.hh:18
std::string m_className
Definition: ViewBase.hh:78
TPad * GetCanvas()
Retrieve the canvas.
Definition: ViewBase.cc:74
void SetRangeY(const double ymin, const double ymax, const bool logscale=false)
Set the range of the function (velocity etc.) to be plotted.
Definition: ViewMedium.cc:87
void SetRangeB(const double bmin, const double bmax, const bool logscale)
Set the limits of the magnetic field.
Definition: ViewMedium.cc:65
void SetMedium(Medium *m)
Set the medium from which to retrieve the transport coefficients.
Definition: ViewMedium.cc:44
void PlotElectronCrossSections()
Definition: ViewMedium.cc:99
void SetRangeE(const double emin, const double emax, const bool logscale)
Set the limits of the electric field.
Definition: ViewMedium.cc:53
void SetRangeA(const double amin, const double amax, const bool logscale)
Set the limits of the angle between electric and magnetic field.
Definition: ViewMedium.cc:76
ViewMedium()
Constructor.
Definition: ViewMedium.cc:42
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384