Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewFEMesh.cc
Go to the documentation of this file.
1#include <algorithm>
2#include <iomanip>
3#include <iostream>
4#include <sstream>
5
6#include <TH1F.h>
7#include <TPolyLine.h>
8#include <TGraph.h>
9
13#include "Garfield/Plotting.hh"
14#include "Garfield/Random.hh"
16
17namespace Garfield {
18
19ViewFEMesh::ViewFEMesh() : ViewBase("ViewFEMesh") {}
20
22 if (!cmp) {
23 std::cerr << m_className << "::SetComponent: Null pointer.\n";
24 return;
25 }
26
27 m_component = cmp;
28}
29
30// The plotting functionality here is ported from Garfield
31// with some inclusion of code from ViewCell.cc
33 if (!m_component) {
34 std::cerr << m_className << "::Plot: Component is not defined.\n";
35 return false;
36 }
37
38 double pmin = 0., pmax = 0.;
39 if (!m_component->GetVoltageRange(pmin, pmax)) {
40 std::cerr << m_className << "::Plot: Component is not ready.\n";
41 return false;
42 }
43
44 if (!GetPlotLimits()) return false;
45
46 auto pad = GetCanvas();
47 pad->cd();
48 if (!RangeSet(pad)) {
50 }
51
52 if (m_drawAxes) {
53 if (!m_xaxis && !m_yaxis) {
54 // Draw default axes.
55 auto frame = pad->DrawFrame(m_xMinPlot, m_yMinPlot,
57 if (m_xaxisTitle.empty()) {
58 frame->GetXaxis()->SetTitle(LabelX().c_str());
59 } else {
60 frame->GetXaxis()->SetTitle(m_xaxisTitle.c_str());
61 }
62 if (m_yaxisTitle.empty()) {
63 frame->GetYaxis()->SetTitle(LabelY().c_str());
64 } else {
65 frame->GetYaxis()->SetTitle(m_yaxisTitle.c_str());
66 }
67 } else {
68 // Draw custom axes.
69 if (m_xaxis) m_xaxis->Draw();
70 if (m_yaxis) m_yaxis->Draw();
71 }
72 }
73
74 // Plot the mesh elements.
75 ComponentCST* componentCST = dynamic_cast<ComponentCST*>(m_component);
76 if (componentCST) {
77 std::cout << m_className << "::Plot: CST component. Calling DrawCST.\n";
78 DrawCST(componentCST);
79 } else {
80 DrawElements();
81 }
82
83 // If we have an associated ViewDrift object, plot the drift lines.
84 if (m_viewDrift) {
85 // Plot a 2D projection of the drift line.
86 for (const auto& dline : m_viewDrift->m_driftLines) {
87 TGraph gr;
88 if (dline.second == ViewDrift::Particle::Electron) {
89 gr.SetLineColor(kOrange - 3);
90 } else {
91 gr.SetLineColor(kRed + 1);
92 }
93 std::vector<float> xgr;
94 std::vector<float> ygr;
95 // Loop over the points.
96 for (const auto& point : dline.first) {
97 // Project this point onto the plane.
98 float xp = 0., yp = 0.;
99 ToPlane(point[0], point[1], point[2], xp, yp);
100 // Add this point if it is within the view.
101 if (InView(xp, yp)) {
102 xgr.push_back(xp);
103 ygr.push_back(yp);
104 }
105 }
106 if (!xgr.empty()) {
107 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "lsame");
108 }
109 }
110 }
111
112 if (m_drawViewRegion && !m_viewRegionX.empty()) {
113 TPolyLine poly;
114 poly.SetLineColor(kSpring + 4);
115 poly.SetLineWidth(3);
116 std::vector<double> xv = m_viewRegionX;
117 std::vector<double> yv = m_viewRegionY;
118 // Close the polygon.
119 xv.push_back(m_viewRegionX[0]);
120 yv.push_back(m_viewRegionY[0]);
121 poly.DrawPolyLine(xv.size(), xv.data(), yv.data(), "same");
122 }
123 gPad->Update();
124 // Draw axes again so they are on top
125 gPad->RedrawAxis("g");
126 return true;
127}
128
129bool ViewFEMesh::GetPlotLimits() {
130
131 if (m_userPlotLimits) {
132 std::vector<double> xp = {m_xMinPlot, m_xMinPlot, m_xMaxPlot, m_xMaxPlot};
133 std::vector<double> yp = {m_yMinPlot, m_yMaxPlot, m_yMaxPlot, m_yMinPlot};
134 std::vector<double> xg(4, 0.);
135 std::vector<double> yg(4, 0.);
136 std::vector<double> zg(4, 0.);
137 for (size_t i = 0; i < 4; ++i) {
138 xg[i] = m_proj[0][0] * xp[i] + m_proj[1][0] * yp[i] + m_proj[2][0];
139 yg[i] = m_proj[0][1] * xp[i] + m_proj[1][1] * yp[i] + m_proj[2][1];
140 zg[i] = m_proj[0][2] * xp[i] + m_proj[1][2] * yp[i] + m_proj[2][2];
141 }
142 m_xMinBox = *std::min_element(xg.begin(), xg.end());
143 m_xMaxBox = *std::max_element(xg.begin(), xg.end());
144 m_yMinBox = *std::min_element(yg.begin(), yg.end());
145 m_yMaxBox = *std::max_element(yg.begin(), yg.end());
146 m_zMinBox = *std::min_element(zg.begin(), zg.end());
147 m_zMaxBox = *std::max_element(zg.begin(), zg.end());
148 m_viewRegionX = xp;
149 m_viewRegionY = yp;
150 return true;
151 }
152 if (!m_userBox) {
153 // If not set by the user, get the bounding box of the component.
154 if (!m_component) return false;
155 if (!m_component->GetBoundingBox(m_xMinBox, m_yMinBox, m_zMinBox,
157 std::cerr << m_className << "::GetPlotLimits:\n"
158 << " Bounding box of the component is not defined.\n"
159 << " Please set the limits explicitly (SetArea).\n";
160 return false;
161 }
162 if (std::isinf(m_xMinBox) || std::isinf(m_xMaxBox) ||
163 std::isinf(m_yMinBox) || std::isinf(m_yMaxBox) ||
164 std::isinf(m_zMinBox) || std::isinf(m_zMaxBox)) {
165 double x0 = 0., y0 = 0., z0 = 0.;
166 double x1 = 0., y1 = 0., z1 = 0.;
167 if (!m_component->GetElementaryCell(x0, y0, z0, x1, y1, z1)) {
168 std::cerr << m_className << "::GetPlotLimits:\n"
169 << " Cell boundaries are not defined.\n"
170 << " Please set the limits explicitly (SetArea).\n";
171 }
172 if (std::isinf(m_xMinBox) || std::isinf(m_xMaxBox)) {
173 m_xMinBox = x0;
174 m_xMaxBox = x1;
175 }
176 if (std::isinf(m_yMinBox) || std::isinf(m_yMaxBox)) {
177 m_yMinBox = y0;
178 m_yMaxBox = y1;
179 }
180 if (std::isinf(m_zMinBox) || std::isinf(m_zMaxBox)) {
181 m_zMinBox = z0;
182 m_zMaxBox = z1;
183 }
184 }
185 }
186 // Determine the intersection of the viewing plane and the box.
187 double xmin = 0., xmax = 0., ymin = 0., ymax = 0.;
188 IntersectPlaneArea(xmin, ymin, xmax, ymax);
189 if (m_viewRegionX.empty()) {
190 std::cerr << m_className << "::GetPlotLimits: Empty view.\n"
191 << " Make sure the viewing plane (SetPlane)\n"
192 << " intersects with the bounding box.\n";
193 return false;
194 }
195 m_xMinPlot = xmin;
196 m_xMaxPlot = xmax;
197 m_yMinPlot = ymin;
198 m_yMaxPlot = ymax;
199 return true;
200}
201
202void ViewFEMesh::SetPlane(const double fx, const double fy, const double fz,
203 const double x0, const double y0, const double z0) {
204 if (fy * fy + fz * fz > 0) {
205 SetPlane(fx, fy, fz, x0, y0, z0, 1, 0, 0);
206 } else {
207 SetPlane(fx, fy, fz, x0, y0, z0, 0, 1, 0);
208 }
209}
210
211void ViewFEMesh::SetPlane(const double fx, const double fy, const double fz,
212 const double x0, const double y0, const double z0,
213 const double hx, const double hy, const double hz) {
214 ViewBase::SetPlane(fx, fy, fz, x0, y0, z0, hx, hy, hz);
215}
216
217// Set the x-axis.
218void ViewFEMesh::SetXaxis(TGaxis* ax) { m_xaxis = ax; }
219
220// Set the y-axis.
221void ViewFEMesh::SetYaxis(TGaxis* ay) { m_yaxis = ay; }
222
223// Create default axes
225 // Create a new x and y axis.
226 if (!GetPlotLimits()) {
227 std::cerr << m_className << "::CreateDefaultAxes:\n"
228 << " Cannot determine the axis limits.\n";
229 return;
230 }
231 const double dx = std::abs(m_xMaxPlot - m_xMinPlot) * 0.1;
232 const double dy = std::abs(m_yMaxPlot - m_yMinPlot) * 0.1;
233 const double x0 = m_xMinPlot + dx;
234 const double y0 = m_yMinPlot + dy;
235 const double x1 = m_xMaxPlot - dx;
236 const double y1 = m_yMaxPlot - dy;
237 m_xaxis = new TGaxis(x0, y0, x1, y0, x0, x1, 2405, "x");
238 m_yaxis = new TGaxis(x0, y0, x0, y1, y0, y1, 2405, "y");
239
240 // Label sizes
241 m_xaxis->SetLabelSize(0.025);
242 m_yaxis->SetLabelSize(0.025);
243
244 // Titles
245 m_xaxis->SetTitleSize(0.03);
246 m_xaxis->SetTitle(LabelX().c_str());
247 m_yaxis->SetTitleSize(0.03);
248 m_yaxis->SetTitle(LabelY().c_str());
249}
250
251// Use ROOT plotting functions to draw the mesh elements on the canvas.
252// General methodology ported from Garfield
253void ViewFEMesh::DrawElements() {
254 // Get the map boundaries from the component.
255 double mapxmax = m_component->m_mapmax[0];
256 double mapxmin = m_component->m_mapmin[0];
257 double mapymax = m_component->m_mapmax[1];
258 double mapymin = m_component->m_mapmin[1];
259 double mapzmax = m_component->m_mapmax[2];
260 double mapzmin = m_component->m_mapmin[2];
261
262 // Get the periodicities.
263 double sx = mapxmax - mapxmin;
264 double sy = mapymax - mapymin;
265 double sz = mapzmax - mapzmin;
266 const bool perX =
267 m_component->m_periodic[0] || m_component->m_mirrorPeriodic[0];
268 const bool perY =
269 m_component->m_periodic[1] || m_component->m_mirrorPeriodic[1];
270 const bool perZ =
271 m_component->m_periodic[2] || m_component->m_mirrorPeriodic[2];
272
273 // Get the plane information.
274 const double fx = m_plane[0];
275 const double fy = m_plane[1];
276 const double fz = m_plane[2];
277 const double dist = m_plane[3];
278
279 // Construct single-column matrix for use as coordinate vector.
280 TMatrixD xMat(3, 1);
281
282 // Determine the number of periods present in the cell.
283 const int nMinX = perX ? int(m_xMinBox / sx) - 1 : 0;
284 const int nMaxX = perX ? int(m_xMaxBox / sx) + 1 : 0;
285 const int nMinY = perY ? int(m_yMinBox / sy) - 1 : 0;
286 const int nMaxY = perY ? int(m_yMaxBox / sy) + 1 : 0;
287 const int nMinZ = perZ ? int(m_zMinBox / sz) - 1 : 0;
288 const int nMaxZ = perZ ? int(m_zMaxBox / sz) + 1 : 0;
289
290 bool cst = false;
291 if (dynamic_cast<ComponentCST*>(m_component)) cst = true;
292
293 // Loop over all elements.
294 const auto nElements = m_component->GetNumberOfElements();
295 for (size_t i = 0; i < nElements; ++i) {
296 size_t mat = 0;
297 bool driftmedium = false;
298 std::vector<size_t> nodes;
299 if (!m_component->GetElement(i, mat, driftmedium, nodes)) continue;
300 // Do not plot the drift medium.
301 if (driftmedium && !m_plotMeshBorders) continue;
302 // Do not create polygons for disabled materials.
303 if (m_disabledMaterial[mat]) continue;
304
305 TGraph gr;
306 const short col = m_colorMap.count(mat) != 0 ? m_colorMap[mat] : 1;
307 gr.SetLineColor(col);
308 if (m_colorMap_fill.count(mat) != 0) {
309 gr.SetFillColor(m_colorMap_fill[mat]);
310 } else {
311 gr.SetFillColor(col);
312 }
313 gr.SetLineWidth(3);
314 std::string opt = "";
315 if (m_plotMeshBorders || !m_fillMesh) opt += "l";
316 if (m_fillMesh) opt += "f";
317 opt += "same";
318
319 // Get the vertex coordinates in the basic cell.
320 std::vector<double> vx0;
321 std::vector<double> vy0;
322 std::vector<double> vz0;
323 const size_t nNodes = nodes.size();
324 for (size_t j = 0; j < nNodes; ++j) {
325 double xn = 0., yn = 0., zn = 0.;
326 if (!m_component->GetNode(nodes[j], xn, yn, zn)) continue;
327 vx0.push_back(xn);
328 vy0.push_back(yn);
329 vz0.push_back(zn);
330 }
331 if (vx0.size() != nNodes) {
332 std::cerr << m_className << "::DrawElements:\n"
333 << " Error retrieving nodes of element " << i << ".\n";
334 continue;
335 }
336 // Coordinates of vertices
337 std::vector<double> vx(nNodes, 0.);
338 std::vector<double> vy(nNodes, 0.);
339 std::vector<double> vz(nNodes, 0.);
340 // Loop over the periodicities in x.
341 for (int nx = nMinX; nx <= nMaxX; nx++) {
342 const double dx = sx * nx;
343 // Determine the x-coordinates of the vertices.
344 if (m_component->m_mirrorPeriodic[0] && nx != 2 * (nx / 2)) {
345 for (size_t j = 0; j < nNodes; ++j) {
346 vx[j] = mapxmin + (mapxmax - vx0[j]) + dx;
347 }
348 } else {
349 for (size_t j = 0; j < nNodes; ++j) {
350 vx[j] = vx0[j] + dx;
351 }
352 }
353
354 // Loop over the periodicities in y.
355 for (int ny = nMinY; ny <= nMaxY; ny++) {
356 const double dy = sy * ny;
357 // Determine the y-coordinates of the vertices.
358 if (m_component->m_mirrorPeriodic[1] && ny != 2 * (ny / 2)) {
359 for (size_t j = 0; j < nNodes; ++j) {
360 vy[j] = mapymin + (mapymax - vy0[j]) + dy;
361 }
362 } else {
363 for (size_t j = 0; j < nNodes; ++j) {
364 vy[j] = vy0[j] + dy;
365 }
366 }
367
368 // Loop over the periodicities in z.
369 for (int nz = nMinZ; nz <= nMaxZ; nz++) {
370 const double dz = sz * nz;
371 // Determine the z-coordinates of the vertices.
372 if (m_component->m_mirrorPeriodic[2] && nz != 2 * (nz / 2)) {
373 for (size_t j = 0; j < nNodes; ++j) {
374 vz[j] = mapzmin + (mapzmax - vz0[j]) + dz;
375 }
376 } else {
377 for (size_t j = 0; j < nNodes; ++j) {
378 vz[j] = vz0[j] + dz;
379 }
380 }
381
382 // Store the x and y coordinates of the relevant mesh vertices.
383 std::vector<double> vX;
384 std::vector<double> vY;
385
386 // Value used to determine whether a vertex is in the plane.
387 const double pcf = std::max(
388 {std::abs(vx[0]), std::abs(vy[0]), std::abs(vz[0]),
389 std::abs(fx), std::abs(fy), std::abs(fz), std::abs(dist)});
390 const double tol = 1.e-4 * pcf;
391 // First isolate the vertices that are in the viewing plane.
392 std::vector<bool> in(nNodes, false);
393 int cnt = 0;
394 for (size_t j = 0; j < nNodes; ++j) {
395 const double d = fx * vx[j] + fy * vy[j] + fz * vz[j] - dist;
396 if (std::abs(d) < tol) {
397 // Point is in the plane.
398 in[j] = true;
399 // Calculate the planar coordinates.
400 double xp = 0., yp = 0.;
401 ToPlane(vx[j], vy[j], vz[j], xp, yp);
402 vX.push_back(xp);
403 vY.push_back(yp);
404 } else {
405 if (d > 0.) {
406 cnt += 1;
407 } else {
408 cnt -= 1;
409 }
410 }
411 }
412 // Stop if all points are on the same side of the plane.
413 if (std::abs(cnt) == (int)nNodes) continue;
414 // Cut the sides that are not in the plane.
415 if (cst) {
416 const std::array<std::array<unsigned int, 3>, 8> neighbours = {{
417 {1, 2, 4}, {0, 3, 5}, {0, 3, 6}, {1, 2, 7},
418 {0, 5, 6}, {1, 4, 7}, {2, 4, 7}, {3, 5, 6}
419 }};
420 for (size_t j = 0; j < nNodes; ++j) {
421 for (unsigned int k : neighbours[j]) {
422 if (in[j] || in[k]) continue;
423 if (PlaneCut(vx[j], vy[j], vz[j],
424 vx[k], vy[k], vz[k], xMat)) {
425 vX.push_back(xMat(0, 0));
426 vY.push_back(xMat(1, 0));
427 }
428 }
429 }
430 } else {
431 // Tetrahedron.
432 for (size_t j = 0; j < nNodes; ++j) {
433 for (size_t k = j + 1; k < nNodes; ++k) {
434 if (in[j] || in[k]) continue;
435 if (PlaneCut(vx[j], vy[j], vz[j],
436 vx[k], vy[k], vz[k], xMat)) {
437 vX.push_back(xMat(0, 0));
438 vY.push_back(xMat(1, 0));
439 }
440 }
441 }
442 }
443 if (vX.size() < 3) continue;
444
445 // Eliminate crossings of the polygon lines
446 // (known as "butterflies" in Garfield).
447 RemoveCrossings(vX, vY);
448
449 // Create vectors to store the clipped polygon.
450 std::vector<double> cX;
451 std::vector<double> cY;
452
453 // Clip the polygon to the view area.
454 ClipToView(vX, vY, cX, cY);
455
456 // If we have more than 2 vertices, add the polygon.
457 if (cX.size() <= 2) continue;
458 // Again eliminate crossings of the polygon lines.
459 RemoveCrossings(cX, cY);
460
461 // Draw the polygon.
462 std::vector<float> xgr(cX.begin(), cX.end());
463 std::vector<float> ygr(cY.begin(), cY.end());
464 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), opt.c_str());
465 } // end z-periodicity loop
466 } // end y-periodicity loop
467 } // end x-periodicity loop
468 } // end loop over elements
469
470}
471
472void ViewFEMesh::DrawCST(ComponentCST* cst) {
473 /*The method is based on ViewFEMesh::Draw, thus the first part is copied from
474 * there.
475 * At the moment only x-y, x-z, and y-z are available due to the simple
476 * implementation.
477 * The advantage of this method is that there is no element loop and thus it
478 * is much
479 * faster.
480 */
481
482 // Helper struct.
483 struct PolygonInfo {
484 double p1[2];
485 double p2[2];
486 double p3[2];
487 double p4[2];
488 int element;
489 };
490
491 // Get the map boundaries from the component
492 double mapxmax = m_component->m_mapmax[0];
493 double mapxmin = m_component->m_mapmin[0];
494 double mapymax = m_component->m_mapmax[1];
495 double mapymin = m_component->m_mapmin[1];
496 double mapzmax = m_component->m_mapmax[2];
497 double mapzmin = m_component->m_mapmin[2];
498
499 // Get the periodicities.
500 double sx = mapxmax - mapxmin;
501 double sy = mapymax - mapymin;
502 double sz = mapzmax - mapzmin;
503 const bool perX =
504 m_component->m_periodic[0] || m_component->m_mirrorPeriodic[0];
505 const bool perY =
506 m_component->m_periodic[1] || m_component->m_mirrorPeriodic[1];
507 const bool perZ =
508 m_component->m_periodic[2] || m_component->m_mirrorPeriodic[2];
509
510 // Determine the number of periods present in the cell.
511 const int nMinX = perX ? int(m_xMinBox / sx) - 1 : 0;
512 const int nMaxX = perX ? int(m_xMaxBox / sx) + 1 : 0;
513 const int nMinY = perY ? int(m_yMinBox / sy) - 1 : 0;
514 const int nMaxY = perY ? int(m_yMaxBox / sy) + 1 : 0;
515 const int nMinZ = perZ ? int(m_zMinBox / sz) - 1 : 0;
516 const int nMaxZ = perZ ? int(m_zMaxBox / sz) + 1 : 0;
517
518 std::vector<PolygonInfo> elements;
519 int nMinU = 0, nMaxU = 0, nMinV = 0, nMaxV = 0;
520 double mapumin = 0., mapumax = 0., mapvmin = 0., mapvmax = 0.;
521 double su = 0., sv = 0.;
522 bool mirroru = false, mirrorv = false;
523 double uMin, vMin, uMax, vMax;
524 unsigned int n_x, n_y, n_z;
525 cst->GetNumberOfMeshLines(n_x, n_y, n_z);
526
527 const double fx = m_plane[0];
528 const double fy = m_plane[1];
529 const double fz = m_plane[2];
530 if (fx == 0 && fy == 0 && fz == 1) {
531 // xy view
532 std::cout << m_className << "::DrawCST: Creating x-y mesh view.\n";
533 // Calculate the z position.
534 unsigned int i, j, z;
535 const double z0 = m_plane[3] * fz;
536 if (!cst->Coordinate2Index(0, 0, z0, i, j, z)) {
537 std::cerr << " Could not determine the z-index of the plane.\n";
538 return;
539 }
540 std::cout << " The z-index of the plane is " << z << ".\n";
541 nMinU = nMinX;
542 nMaxU = nMaxX;
543 nMinV = nMinY;
544 nMaxV = nMaxY;
545 uMin = m_xMinBox;
546 uMax = m_xMaxBox;
547 vMin = m_yMinBox;
548 vMax = m_yMaxBox;
549
550 mapumin = mapxmin;
551 mapumax = mapxmax;
552 mapvmin = mapymin;
553 mapvmax = mapymax;
554 su = sx;
555 sv = sy;
556 mirroru = perX;
557 mirrorv = perY;
558 for (unsigned int y = 0; y < (n_y - 1); y++) {
559 for (unsigned int x = 0; x < (n_x - 1); x++) {
560 auto elem = cst->Index2Element(x, y, z);
561 double e_xmin, e_xmax, e_ymin, e_ymax, e_zmin, e_zmax;
562 cst->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
563 e_zmin, e_zmax);
564 PolygonInfo tmp_info;
565 tmp_info.element = elem;
566 tmp_info.p1[0] = e_xmin;
567 tmp_info.p2[0] = e_xmax;
568 tmp_info.p3[0] = e_xmax;
569 tmp_info.p4[0] = e_xmin;
570 tmp_info.p1[1] = e_ymin;
571 tmp_info.p2[1] = e_ymin;
572 tmp_info.p3[1] = e_ymax;
573 tmp_info.p4[1] = e_ymax;
574 elements.push_back(std::move(tmp_info));
575 }
576 }
577 } else if (fx == 0 && fy == -1 && fz == 0) {
578 // xz-view
579 std::cout << m_className << "::DrawCST: Creating x-z mesh view.\n";
580 // Calculate the y position.
581 unsigned int i = 0, j = 0, y = 0;
582 const double y0 = m_plane[3] * fy;
583 if (!cst->Coordinate2Index(0, y0, 0, i, y, j)) {
584 std::cerr << " Could not determine the y-index of the plane.\n";
585 return;
586 }
587 std::cout << " The y-index of the plane is " << y << ".\n";
588
589 nMinU = nMinX;
590 nMaxU = nMaxX;
591 nMinV = nMinZ;
592 nMaxV = nMaxZ;
593 uMin = m_xMinBox;
594 uMax = m_xMaxBox;
595 vMin = m_zMinBox;
596 vMax = m_zMaxBox;
597
598 mapumin = mapxmin;
599 mapumax = mapxmax;
600 mapvmin = mapzmin;
601 mapvmax = mapzmax;
602 su = sx;
603 sv = sz;
604 mirroru = perX;
605 mirrorv = perZ;
606 for (unsigned int z = 0; z < (n_z - 1); z++) {
607 for (unsigned int x = 0; x < (n_x - 1); x++) {
608 auto elem = cst->Index2Element(x, y, z);
609 double e_xmin, e_xmax, e_ymin, e_ymax, e_zmin, e_zmax;
610 cst->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
611 e_zmin, e_zmax);
612 PolygonInfo tmp_info;
613 tmp_info.element = elem;
614 tmp_info.p1[0] = e_xmin;
615 tmp_info.p2[0] = e_xmax;
616 tmp_info.p3[0] = e_xmax;
617 tmp_info.p4[0] = e_xmin;
618 tmp_info.p1[1] = e_zmin;
619 tmp_info.p2[1] = e_zmin;
620 tmp_info.p3[1] = e_zmax;
621 tmp_info.p4[1] = e_zmax;
622 elements.push_back(std::move(tmp_info));
623 }
624 }
625 } else if (fx == -1 && fy == 0 && fz == 0) {
626 // yz-view
627 std::cout << m_className << "::DrawCST: Creating z-y mesh view.\n";
628 // Calculate the x position.
629 unsigned int i, j, x;
630 const double x0 = m_plane[3] * fx;
631 if (!cst->Coordinate2Index(x0, 0, 0, x, i, j)) {
632 std::cerr << " Could not determine the x-index of the plane.\n";
633 return;
634 }
635 std::cout << " The x-index of the plane is " << x << ".\n";
636 nMinU = nMinZ;
637 nMaxU = nMaxZ;
638 nMinV = nMinY;
639 nMaxV = nMaxY;
640 uMin = m_yMinBox;
641 uMax = m_yMaxBox;
642 vMin = m_zMinBox;
643 vMax = m_zMaxBox;
644
645 mapumin = mapzmin;
646 mapumax = mapzmax;
647 mapvmin = mapymin;
648 mapvmax = mapymax;
649 su = sz;
650 sv = sy;
651 mirroru = perZ;
652 mirrorv = perY;
653 for (unsigned int z = 0; z < (n_z - 1); z++) {
654 for (unsigned int y = 0; y < (n_y - 1); y++) {
655 auto elem = cst->Index2Element(x, y, z);
656 double e_xmin, e_xmax, e_ymin, e_ymax, e_zmin, e_zmax;
657 cst->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
658 e_zmin, e_zmax);
659 PolygonInfo tmp_info;
660 tmp_info.element = elem;
661 tmp_info.p1[0] = e_zmin;
662 tmp_info.p2[0] = e_zmax;
663 tmp_info.p3[0] = e_zmax;
664 tmp_info.p4[0] = e_zmin;
665 tmp_info.p1[1] = e_ymin;
666 tmp_info.p2[1] = e_ymin;
667 tmp_info.p3[1] = e_ymax;
668 tmp_info.p4[1] = e_ymax;
669 // Add the polygon to the mesh
670 elements.push_back(std::move(tmp_info));
671 }
672 }
673 } else {
674 std::cerr << m_className << "::DrawCST:\n";
675 std::cerr << " The given plane name is not known.\n";
676 std::cerr << " Please choose one of the following: xy, xz, yz.\n";
677 return;
678 }
679
680 std::cout << m_className << "::DrawCST:\n " << elements.size()
681 << " elements in the projection of the unit cell.\n";
682 for (const auto& element : elements) {
683 size_t mat = 0;
684 bool driftmedium = false;
685 std::vector<size_t> nodes;
686 cst->GetElement(element.element, mat, driftmedium, nodes);
687 // Do not plot the drift medium.
688 if (driftmedium && !m_plotMeshBorders) continue;
689 // Do not create polygons for disabled materials.
690 if (m_disabledMaterial[mat]) continue;
691 TGraph gr;
692 const short col = m_colorMap.count(mat) > 0 ? m_colorMap[mat] : 1;
693 gr.SetLineColor(col);
694 if (m_colorMap_fill.count(mat) > 0) {
695 gr.SetFillColor(m_colorMap_fill[mat]);
696 } else {
697 gr.SetFillColor(col);
698 }
699 if (m_plotMeshBorders)
700 gr.SetLineWidth(3);
701 else
702 gr.SetLineWidth(1);
703
704 for (int nu = nMinU; nu <= nMaxU; nu++) {
705 for (int nv = nMinV; nv <= nMaxV; nv++) {
706 // Add 4 points of the square
707 float tmp_u[4], tmp_v[4];
708 if (mirroru && nu != 2 * (nu / 2)) {
709 // nu is odd
710 tmp_u[0] = mapumin + (mapumax - element.p1[0]) + su * nu;
711 tmp_u[1] = mapumin + (mapumax - element.p2[0]) + su * nu;
712 tmp_u[2] = mapumin + (mapumax - element.p3[0]) + su * nu;
713 tmp_u[3] = mapumin + (mapumax - element.p4[0]) + su * nu;
714 } else {
715 // nu is even
716 tmp_u[0] = element.p1[0] + su * nu;
717 tmp_u[1] = element.p2[0] + su * nu;
718 tmp_u[2] = element.p3[0] + su * nu;
719 tmp_u[3] = element.p4[0] + su * nu;
720 }
721 if (mirrorv && nv != 2 * (nv / 2)) {
722 tmp_v[0] = mapvmin + (mapvmax - element.p1[1]) + sv * nv;
723 tmp_v[1] = mapvmin + (mapvmax - element.p2[1]) + sv * nv;
724 tmp_v[2] = mapvmin + (mapvmax - element.p3[1]) + sv * nv;
725 tmp_v[3] = mapvmin + (mapvmax - element.p4[1]) + sv * nv;
726 } else {
727 tmp_v[0] = element.p1[1] + sv * nv;
728 tmp_v[1] = element.p2[1] + sv * nv;
729 tmp_v[2] = element.p3[1] + sv * nv;
730 tmp_v[3] = element.p4[1] + sv * nv;
731 }
732 if (tmp_u[0] < uMin || tmp_u[1] > uMax || tmp_v[0] < vMin ||
733 tmp_v[2] > vMax) {
734 continue;
735 }
736 std::string opt = "";
737 if (m_plotMeshBorders || !m_fillMesh) opt += "l";
738 if (m_fillMesh) opt += "f";
739 opt += "same";
740 gr.DrawGraph(4, tmp_u, tmp_v, opt.c_str());
741 }
742 }
743 }
744
745}
746
747// Removes duplicate points and line crossings by correctly ordering
748// the points in the provided vectors.
749//
750// NOTE: This is a 2D version of the BUTFLD method in Garfield. It
751// follows the same general algorithm.
752//
753// TODO: there is an algorithm which always sorts points correctly using cross
754// product, see IntersectPlaneArea.
755void ViewFEMesh::RemoveCrossings(std::vector<double>& x,
756 std::vector<double>& y) {
757 // Determine element dimensions
758 double xmin = x[0], xmax = x[0];
759 double ymin = y[0], ymax = y[0];
760 for (int i = 1; i < (int)x.size(); i++) {
761 if (x[i] < xmin) xmin = x[i];
762 if (x[i] > xmax) xmax = x[i];
763 if (y[i] < ymin) ymin = y[i];
764 if (y[i] > ymax) ymax = y[i];
765 }
766
767 // First remove duplicate points
768 double xtol = 1e-10 * std::abs(xmax - xmin);
769 double ytol = 1e-10 * std::abs(ymax - ymin);
770 for (int i = 0; i < (int)x.size(); i++) {
771 for (int j = i + 1; j < (int)x.size(); j++) {
772 if (std::abs(x[i] - x[j]) < xtol && std::abs(y[i] - y[j]) < ytol) {
773 x.erase(x.begin() + j);
774 y.erase(y.begin() + j);
775 j--;
776 }
777 }
778 }
779
780 // No crossings with 3 points or less
781 if (x.size() <= 3) return;
782
783 // Save the polygon size so it is easily accessible
784 int NN = x.size();
785
786 // Keep track of the number of attempts
787 int attempts = 0;
788
789 // Exchange points until crossings are eliminated or we have attempted NN
790 // times
791 bool crossings = true;
792 while (crossings && (attempts < NN)) {
793 // Assume we are done after this attempt.
794 crossings = false;
795
796 for (int i = 1; i <= NN; i++) {
797 for (int j = i + 2; j <= NN; j++) {
798 // End the j-loop if we have surpassed N and wrapped around to i
799 if ((j + 1) > NN && 1 + (j % NN) >= i) break;
800 // Otherwise, detect crossings and attempt to eliminate them.
801 // Determine if we have a crossing.
802 double xc = 0., yc = 0.;
803 if (!LinesCrossed(x[(i - 1) % NN], y[(i - 1) % NN], x[i % NN],
804 y[i % NN], x[(j - 1) % NN], y[(j - 1) % NN],
805 x[j % NN], y[j % NN], xc, yc)) {
806 continue;
807 }
808 // Swap each point from i towards j with each corresponding point
809 // from j towards i.
810 for (int k = 1; k <= (j - i) / 2; k++) {
811 double xs = x[(i + k - 1) % NN];
812 double ys = y[(i + k - 1) % NN];
813 x[(i + k - 1) % NN] = x[(j - k) % NN];
814 y[(i + k - 1) % NN] = y[(j - k) % NN];
815 x[(j - k) % NN] = xs;
816 y[(j - k) % NN] = ys;
817
818 // Force another attempt
819 crossings = true;
820 }
821 } // end loop over j
822 } // end loop over i
823
824 // Increment the number of attempts
825 attempts++;
826
827 } // end while(crossings)
828
829 if (attempts > NN) {
830 std::cerr << m_className << "::RemoveCrossings:\n Warning: "
831 << "Maximum attempts reached. Crossings not removed.\n";
832 }
833}
834
835/// Return true if the specified point is in the view region.
836bool ViewFEMesh::InView(const double x, const double y) const {
837 // Test whether this vertex is inside the view.
838 if (m_userPlotLimits) {
839 return (x >= m_xMinPlot && x <= m_xMaxPlot &&
840 y >= m_yMinPlot && y <= m_yMaxPlot);
841 }
842 bool edge = false;
843 return IsInPolygon(x, y, m_viewRegionX, m_viewRegionY, edge);
844}
845
846//
847// Determines whether the line connecting points (x1,y1) and (x2,y2)
848// and the line connecting points (u1,v1) and (u2,v2) cross somewhere
849// between the 4 points. Sets the crossing point in (xc, yc).
850//
851// Ported from Garfield function CROSSD
852//
853bool ViewFEMesh::LinesCrossed(double x1, double y1, double x2, double y2,
854 double u1, double v1, double u2, double v2,
855 double& xc, double& yc) const {
856 // Set the tolerances.
857 double xtol = 1.0e-10 * std::max({std::abs(x1), std::abs(x2), std::abs(u1),
858 std::abs(u2)});
859 double ytol = 1.0e-10 * std::max({std::abs(y1), std::abs(y2), std::abs(v1),
860 std::abs(v2)});
861 if (xtol <= 0) xtol = 1.0e-10;
862 if (ytol <= 0) ytol = 1.0e-10;
863
864 // Compute the distances and determinant (dx,dy) x (du,dv).
865 double dy = y2 - y1;
866 double dv = v2 - v1;
867 double dx = x1 - x2;
868 double du = u1 - u2;
869 double det = dy * du - dx * dv;
870
871 // Check for crossing because one of the endpoints is on both lines.
872 if (OnLine(x1, y1, x2, y2, u1, v1)) {
873 xc = u1;
874 yc = v1;
875 return true;
876 } else if (OnLine(x1, y1, x2, y2, u2, v2)) {
877 xc = u2;
878 yc = v2;
879 return true;
880 } else if (OnLine(u1, v1, u2, v2, x1, y1)) {
881 xc = x1;
882 yc = y1;
883 return true;
884 } else if (OnLine(u1, v1, u2, v2, x2, y2)) {
885 xc = x2;
886 yc = y2;
887 return true;
888 }
889 // Check if the lines are parallel (zero determinant).
890 if (std::abs(det) < xtol * ytol) return false;
891 // No special case: compute point of intersection.
892
893 // Solve crossing equations.
894 xc = (du * (x1 * y2 - x2 * y1) - dx * (u1 * v2 - u2 * v1)) / det;
895 yc = ((-1 * dv) * (x1 * y2 - x2 * y1) + dy * (u1 * v2 - u2 * v1)) / det;
896
897 // Determine if this point is on both lines.
898 if (OnLine(x1, y1, x2, y2, xc, yc) && OnLine(u1, v1, u2, v2, xc, yc))
899 return true;
900
901 // The lines do not cross if we have reached this point.
902 return false;
903}
904
905// Determines whether the point (u,v) lies on the line connecting
906// points (x1,y1) and (x2,y2).
907//
908// Ported from Garfield function ONLIND
909//
910bool ViewFEMesh::OnLine(double x1, double y1, double x2, double y2, double u,
911 double v) const {
912 // Set the tolerances
913 double xtol = 1.e-10 * std::max({std::abs(x1), std::abs(x2), std::abs(u)});
914 double ytol = 1.e-10 * std::max({std::abs(y1), std::abs(y2), std::abs(v)});
915 if (xtol <= 0) xtol = 1.0e-10;
916 if (ytol <= 0) ytol = 1.0e-10;
917
918 // To store the coordinates of the comparison point
919 double xc = 0, yc = 0;
920
921 // Check if (u,v) is the point (x1,y1) or (x2,y2)
922 if ((std::abs(x1 - u) <= xtol && std::abs(y1 - v) <= ytol) ||
923 (std::abs(x2 - u) <= xtol && std::abs(y2 - v) <= ytol)) {
924 return true;
925 }
926 // Check if the line is actually a point
927 if (std::abs(x1 - x2) <= xtol && std::abs(y1 - y2) <= ytol) {
928 return false;
929 }
930 // Choose (x1,y1) as starting point if closer to (u,v)
931 if (std::abs(u - x1) + std::abs(v - y1) <
932 std::abs(u - x2) + std::abs(v - y2)) {
933 // Compute the component of the line from (x1,y1) to (u,v)
934 // along the line from (x1,y1) to (x2,y2)
935 double dpar = ((u - x1) * (x2 - x1) + (v - y1) * (y2 - y1)) /
936 ((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1));
937
938 // Determine the point on the line to which to compare (u,v)
939 if (dpar < 0.0) {
940 xc = x1;
941 yc = y1;
942 } else if (dpar > 1.0) {
943 xc = x2;
944 yc = y2;
945 } else {
946 xc = x1 + dpar * (x2 - x1);
947 yc = y1 + dpar * (y2 - y1);
948 }
949 } else {
950 // Choose (x2,y2) as starting point if closer to (u,v)
951 // Compute the component of the line from (x2,y2) to (u,v)
952 // along the line from (x2,y2) to (x1,y1)
953 double dpar = ((u - x2) * (x1 - x2) + (v - y2) * (y1 - y2)) /
954 ((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1));
955
956 // Determine the point on the line to which to compare (u,v)
957 if (dpar < 0.0) {
958 xc = x2;
959 yc = y2;
960 } else if (dpar > 1.0) {
961 xc = x1;
962 yc = y1;
963 } else {
964 xc = x2 + dpar * (x1 - x2);
965 yc = y2 + dpar * (y1 - y2);
966 }
967 }
968
969 // Compare the calculated point to (u,v)
970 if (std::abs(u - xc) < xtol && std::abs(v - yc) < ytol) return true;
971
972 return false;
973}
974
975// Ported from Garfield: determines the point of intersection, in planar
976// coordinates, of a plane with the line connecting multiple points
977// x1,y1,z1;x2,y2,z2: the world coordinates of the two points
978// projMat;planeMat: the projection and plane matrices
979// xMat: the resulting planar coordinates of the intersection point
980bool ViewFEMesh::PlaneCut(double x1, double y1, double z1, double x2, double y2,
981 double z2, TMatrixD& xMat) {
982 // Set up the matrix for cutting edges not in the plane
983 TArrayD dataCut(9);
984 TMatrixD cutMat(3, 3);
985 dataCut[0] = m_proj[0][0];
986 dataCut[1] = m_proj[1][0];
987 dataCut[2] = x1 - x2;
988 dataCut[3] = m_proj[0][1];
989 dataCut[4] = m_proj[1][1];
990 dataCut[5] = y1 - y2;
991 dataCut[6] = m_proj[0][2];
992 dataCut[7] = m_proj[1][2];
993 dataCut[8] = z1 - z2;
994 cutMat.SetMatrixArray(dataCut.GetArray());
995
996 // Calculate the determinant of the cut matrix
997 double cutDet =
998 cutMat(0, 0) *
999 (cutMat(1, 1) * cutMat(2, 2) - cutMat(1, 2) * cutMat(2, 1)) -
1000 cutMat(0, 1) *
1001 (cutMat(1, 0) * cutMat(2, 2) - cutMat(1, 2) * cutMat(2, 0)) +
1002 cutMat(0, 2) *
1003 (cutMat(1, 0) * cutMat(2, 1) - cutMat(1, 1) * cutMat(2, 0));
1004
1005 // Do not proceed if the matrix is singular
1006 if (std::abs(cutDet) < 1e-20) return false;
1007
1008 // Set up a coordinate vector (RHS of equation)
1009 TArrayD dataCoords(3);
1010 TMatrixD coordMat(3, 1);
1011 dataCoords[0] = x1 - m_proj[2][0];
1012 dataCoords[1] = y1 - m_proj[2][1];
1013 dataCoords[2] = z1 - m_proj[2][2];
1014 coordMat.SetMatrixArray(dataCoords.GetArray());
1015
1016 // Invert the cut matrix and multiply to get the solution
1017 cutMat.SetTol(1e-20);
1018 cutMat.Invert();
1019 // Do not proceed if the matrix is singular
1020 if (!cutMat.IsValid()) return false;
1021 xMat = cutMat * coordMat;
1022
1023 // Return success if the plane point is between the two vertices
1024 if (xMat(2, 0) < 0 || xMat(2, 0) > 1) return false;
1025 return true;
1026}
1027
1028// Calculates view region and canvas dimensions based on projection plane
1029// and view area
1030bool ViewFEMesh::IntersectPlaneArea(double& xmin, double& ymin,
1031 double& xmax, double& ymax) {
1032 std::vector<TMatrixD> intersect_points;
1033 m_viewRegionX.clear();
1034 m_viewRegionY.clear();
1035 // Loop over box edges
1036 for (int i0 = 0; i0 < 2; ++i0) {
1037 for (int j0 = 0; j0 < 2; ++j0) {
1038 for (int k0 = 0; k0 < 2; ++k0) {
1039 for (int i1 = i0; i1 < 2; ++i1) {
1040 for (int j1 = j0; j1 < 2; ++j1) {
1041 for (int k1 = k0; k1 < 2; ++k1) {
1042 if (i1 - i0 + j1 - j0 + k1 - k0 != 1) continue;
1043 const double x0 = i0 ? m_xMinBox : m_xMaxBox;
1044 const double y0 = j0 ? m_yMinBox : m_yMaxBox;
1045 const double z0 = k0 ? m_zMinBox : m_zMaxBox;
1046 const double x1 = i1 ? m_xMinBox : m_xMaxBox;
1047 const double y1 = j1 ? m_yMinBox : m_yMaxBox;
1048 const double z1 = k1 ? m_zMinBox : m_zMaxBox;
1049 TMatrixD xMat(3, 1);
1050 if (!PlaneCut(x0, y0, z0, x1, y1, z1, xMat)) continue;
1051 if (m_debug) {
1052 std::cout << m_className << "::IntersectPlaneArea:\n"
1053 << " Intersection of plane at (" << xMat(0, 0)
1054 << ", " << xMat(1, 0) << ", " << xMat(2, 0)
1055 << ") with edge\n ("
1056 << x0 << ", " << y0 << ", " << z0 << ")-("
1057 << x1 << ", " << y1 << ", " << z1 << ")\n";
1058 }
1059 // Do not add same points (the case when plane contains an edge)
1060 bool skip = false;
1061 for (auto& p : intersect_points) {
1062 const double dx = xMat(0, 0) - p(0, 0);
1063 const double dy = xMat(1, 0) - p(1, 0);
1064 if (std::hypot(dx, dy) < 1e-10) {
1065 skip = true;
1066 break;
1067 }
1068 }
1069 if (!skip) intersect_points.push_back(xMat);
1070 }
1071 }
1072 }
1073 }
1074 }
1075 }
1076 if (intersect_points.size() < 3) {
1077 std::cerr << m_className << "::IntersectPlaneArea:\n"
1078 << " WARNING: Empty intersection of view plane with area.\n";
1079 return false;
1080 }
1081 TMatrixD offset = intersect_points[0];
1082 xmin = xmax = intersect_points[0](0, 0);
1083 ymin = ymax = intersect_points[0](1, 0);
1084 // Remove crossings by sorting points rotation-wise.
1085 for (auto& p : intersect_points) p -= offset;
1086 std::sort(intersect_points.begin(), intersect_points.end(),
1087 [](const TMatrixD& a, const TMatrixD& b) -> bool {
1088 double cross_z = a(0, 0) * b(1, 0) - a(1, 0) * b(0, 0);
1089 return cross_z < 0;
1090 });
1091 for (auto& p : intersect_points) {
1092 p += offset;
1093 m_viewRegionX.push_back(p(0, 0));
1094 m_viewRegionY.push_back(p(1, 0));
1095 xmin = std::min(p(0, 0), xmin);
1096 ymin = std::min(p(1, 0), ymin);
1097 xmax = std::max(p(0, 0), xmax);
1098 ymax = std::max(p(1, 0), ymax);
1099 }
1100 return true;
1101}
1102
1103// Ported from Garfield (function INTERD):
1104// Returns true if the point (x,y) is inside of the specified polygon.
1105// x: the x-coordinate
1106// y: the y-coordinate
1107// px: the x-vertices of the polygon
1108// py: the y-vertices of the polygon
1109// edge: a variable set to true if the point is located on the polygon edge
1110bool ViewFEMesh::IsInPolygon(double x, double y,
1111 const std::vector<double>& px,
1112 const std::vector<double>& py, bool& edge) const {
1113 // Get the number and coordinates of the polygon vertices.
1114 const size_t pN = px.size();
1115
1116 // Handle the special case of less than 2 vertices.
1117 if (pN < 2) return false;
1118 // Handle the special case of exactly 2 vertices (a line).
1119 if (pN == 2) return OnLine(px[0], py[0], px[1], py[1], x, y);
1120
1121 // Set the minimum and maximum coordinates of all polygon vertices.
1122 double px_min = px[0], py_min = py[0];
1123 double px_max = px[0], py_max = py[0];
1124 for (size_t i = 0; i < pN; i++) {
1125 px_min = std::min(px_min, px[i]);
1126 py_min = std::min(py_min, py[i]);
1127 px_max = std::max(px_max, px[i]);
1128 py_max = std::max(py_max, py[i]);
1129 }
1130
1131 // Set the tolerances
1132 double xtol = 1.0e-10 * std::max(std::abs(px_min), std::abs(px_max));
1133 double ytol = 1.0e-10 * std::max(std::abs(py_min), std::abs(py_max));
1134 if (xtol <= 0) xtol = 1.0e-10;
1135 if (ytol <= 0) ytol = 1.0e-10;
1136
1137 // If we have essentially one x value, check to see if y is in range.
1138 if (std::abs(px_max - px_min) < xtol) {
1139 edge = (y > (py_min - ytol) && y < (py_max + ytol) &&
1140 std::abs(px_max + px_min - 2 * x) < xtol);
1141 return false;
1142 }
1143 // If we have essentially one y value, check to see if x is in range.
1144 if (std::abs(py_max - py_min) < ytol) {
1145 edge = (x > (px_min - xtol) && x < (px_max + xtol) &&
1146 std::abs(py_max + py_min - 2 * y) < ytol);
1147 return false;
1148 }
1149
1150 // Set "infinity" points.
1151 double xinf = px_min - std::abs(px_max - px_min);
1152 double yinf = py_min - std::abs(py_max - py_min);
1153
1154 // Loop until successful or maximum iterations (100) reached.
1155 int niter = 0;
1156 bool done = false;
1157 int ncross = 0;
1158
1159 while (!done && niter < 100) {
1160 // Assume we will finish on this loop.
1161 done = true;
1162
1163 // Loop over all edges, counting the number of edges crossed by a line
1164 // extending from (x, y) to (xinf, yinf).
1165 ncross = 0;
1166 for (size_t i = 0; (done && i < pN); i++) {
1167 // Determine whether the point lies on the edge.
1168 if (OnLine(px[i % pN], py[i % pN], px[(i + 1) % pN], py[(i + 1) % pN], x,
1169 y)) {
1170 edge = true;
1171 return false;
1172 }
1173
1174 // Determine whether this edge is crossed; if so increment the counter.
1175 double xc = 0., yc = 0.;
1176 if (LinesCrossed(x, y, xinf, yinf, px[i % pN], py[i % pN],
1177 px[(i + 1) % pN], py[(i + 1) % pN], xc, yc))
1178 ncross++;
1179
1180 // Ensure this vertex is not crossed by the line from (x,y)
1181 // to (xinf,yinf); if so recompute (xinf,yinf) and start over.
1182 if (OnLine(x, y, xinf, yinf, px[i], py[i])) {
1183 // Recompute (xinf,yinf).
1184 xinf = px_min - RndmUniform() * std::abs(px_max - xinf);
1185 yinf = py_min - RndmUniform() * std::abs(py_max - yinf);
1186
1187 // Start over.
1188 done = false;
1189 niter++;
1190 }
1191 }
1192 }
1193
1194 // If we failed to finish iterating, return false.
1195 if (niter >= 100) {
1196 std::cerr << m_className << "::IsInPolygon: Unable to determine whether ("
1197 << x << ", " << y << ") is inside a polygon. Returning false.\n";
1198 return false;
1199 }
1200
1201 // Point is inside for an odd, nonzero number of crossings.
1202 return (ncross != 2 * (ncross / 2));
1203}
1204
1205// Ported from Garfield (method GRCONV):
1206// Clip the specified polygon to the view region; return the clipped polygon.
1207// px: the x-vertices of the polygon
1208// py: the y-vertices of the polygon
1209// cx: to contain the x-vertices of the clipped polygon
1210// cy: to contain the y-vertices of the clipped polygon
1211void ViewFEMesh::ClipToView(std::vector<double>& px, std::vector<double>& py,
1212 std::vector<double>& cx, std::vector<double>& cy) {
1213 // Get the number and coordinates of the polygon vertices.
1214 int pN = (int)px.size();
1215
1216 // Clear the vectors to contain the final polygon.
1217 cx.clear();
1218 cy.clear();
1219
1220 // Set up the view vertices.
1221 const auto& vx = m_viewRegionX;
1222 const auto& vy = m_viewRegionY;
1223 const int vN = m_viewRegionX.size();
1224
1225 // Do nothing if we have less than 2 points.
1226 if (pN < 2) return;
1227
1228 // Loop over the polygon vertices.
1229 for (int i = 0; i < pN; i++) {
1230 // Flag for skipping check for edge intersection.
1231 bool skip = false;
1232
1233 // Loop over the view vertices.
1234 for (int j = 0; j < vN; j++) {
1235 // Determine whether this vertex lies on a view edge:
1236 // if so add the vertex to the final polygon.
1237 if (OnLine(vx[j % vN], vy[j % vN], vx[(j + 1) % vN], vy[(j + 1) % vN],
1238 px[i], py[i])) {
1239 // Add the vertex.
1240 cx.push_back(px[i]);
1241 cy.push_back(py[i]);
1242
1243 // Skip edge intersection check in this case.
1244 skip = true;
1245 }
1246
1247 // Determine whether a corner of the view area lies on this edge:
1248 // if so add the corner to the final polygon.
1249 if (OnLine(px[i % pN], py[i % pN], px[(i + 1) % pN], py[(i + 1) % pN],
1250 vx[j], vy[j])) {
1251 // Add the vertex.
1252 cx.push_back(vx[j]);
1253 cy.push_back(vy[j]);
1254
1255 // Skip edge intersection check in this case.
1256 skip = true;
1257 }
1258 }
1259
1260 // If we have not skipped the edge intersection check, look for an
1261 // intersection between this edge and the view edges.
1262 if (skip) continue;
1263 // Loop over the view vertices.
1264 for (int j = 0; j < vN; j++) {
1265 // Check for a crossing with this edge;
1266 // if one exists, add the crossing point.
1267 double xc = 0., yc = 0.;
1268 if (LinesCrossed(vx[j % vN], vy[j % vN], vx[(j + 1) % vN],
1269 vy[(j + 1) % vN], px[i % pN], py[i % pN],
1270 px[(i + 1) % pN], py[(i + 1) % pN], xc, yc)) {
1271 // Add a vertex.
1272 cx.push_back(xc);
1273 cy.push_back(yc);
1274 }
1275 }
1276 }
1277
1278 // Find all view field vertices inside the polygon.
1279 for (int j = 0; j < vN; j++) {
1280 // Test whether this vertex is inside the polygon.
1281 // If so, add it to the final polygon.
1282 bool edge = false;
1283 if (IsInPolygon(vx[j], vy[j], px, py, edge)) {
1284 // Add the view vertex.
1285 cx.push_back(vx[j]);
1286 cy.push_back(vy[j]);
1287 }
1288 }
1289
1290 // Find all polygon vertices inside the box.
1291 for (int i = 0; i < pN; i++) {
1292 // Test whether this vertex is inside the view.
1293 // If so, add it to the final polygon.
1294 bool edge = false;
1295 if (IsInPolygon(px[i], py[i], vx, vy, edge)) {
1296 // Add the polygon vertex.
1297 cx.push_back(px[i]);
1298 cy.push_back(py[i]);
1299 }
1300 }
1301}
1302} // namespace Garfield
Base class for components based on finite-element field maps.
virtual bool GetNode(const size_t i, double &x, double &y, double &z) const
bool GetVoltageRange(double &vmin, double &vmax) override
Calculate the voltage range [V].
std::array< double, 3 > m_mapmin
bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override
Get the bounding box coordinates.
virtual size_t GetNumberOfElements() const
Return the number of mesh elements.
bool GetElementaryCell(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override
Get the coordinates of the elementary cell.
bool GetElement(const size_t i, double &vol, double &dmin, double &dmax)
Return the volume and aspect ratio of a mesh element.
std::array< double, 3 > m_mapmax
std::array< bool, 3 > m_mirrorPeriodic
Mirror periodicity in x, y, z.
Definition: Component.hh:346
std::array< bool, 3 > m_periodic
Simple periodicity in x, y, z.
Definition: Component.hh:344
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
bool RangeSet(TPad *)
Definition: ViewBase.cc:83
TPad * GetCanvas()
Retrieve the canvas.
Definition: ViewBase.cc:74
void SetRange(TPad *pad, const double x0, const double y0, const double x1, const double y1)
Definition: ViewBase.cc:93
void ToPlane(const T x, const T y, const T z, T &xp, T &yp) const
Definition: ViewBase.hh:109
virtual void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
Definition: ViewBase.cc:142
void CreateDefaultAxes()
Create a default set of custom-made axes.
Definition: ViewFEMesh.cc:224
ViewFEMesh()
Constructor.
Definition: ViewFEMesh.cc:19
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0) override
Definition: ViewFEMesh.cc:202
bool Plot()
Plot method to be called by user.
Definition: ViewFEMesh.cc:32
void SetComponent(ComponentFieldMap *cmp)
Set the component from which to retrieve the mesh and field.
Definition: ViewFEMesh.cc:21
void SetXaxis(TGaxis *ax)
Definition: ViewFEMesh.cc:218
void SetYaxis(TGaxis *ay)
Definition: ViewFEMesh.cc:221
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14