Garfield++ v2r0
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// Some code was copied/modified from ViewField.cc
2#include <iostream>
3#include <cmath>
4
6#include "ComponentCST.hh"
7#include "Plotting.hh"
8#include "Random.hh"
9#include "ViewFEMesh.hh"
10
11namespace Garfield {
12
14 : m_className("ViewFEMesh"),
15 m_label("Mesh"),
16 m_debug(false),
17 m_fillMesh(false),
18 m_canvas(NULL),
19 m_hasExternalCanvas(false),
20 m_hasUserArea(false),
21 m_xMin(-1.),
22 m_yMin(-1.),
23 m_zMin(-1.),
24 m_xMax(1.),
25 m_yMax(1.),
26 m_zMax(1.),
27 m_component(NULL),
28 m_viewDrift(NULL),
29 m_plotMeshBorders(false),
30 m_xaxis(NULL),
31 m_yaxis(NULL),
32 m_axes(NULL),
33 m_drawAxes(false) {
34
37
38 // Create a blank histogram for the axes.
39 m_axes = new TH2D();
40 m_axes->SetStats(false);
41 m_axes->GetXaxis()->SetTitle("x");
42 m_axes->GetYaxis()->SetTitle("y");
43}
44
46
47 if (!m_hasExternalCanvas && m_canvas) delete m_canvas;
48}
49
51
52 if (!comp) {
53 std::cerr << m_className << "::SetComponent:\n";
54 std::cerr << " Component pointer is null.\n";
55 return;
56 }
57
58 m_component = comp;
59}
60
61void ViewFEMesh::SetCanvas(TCanvas* c) {
62
63 if (!c) return;
64 if (!m_hasExternalCanvas && m_canvas) {
65 delete m_canvas;
66 m_canvas = NULL;
67 }
68 m_canvas = c;
69 m_hasExternalCanvas = true;
70}
71
72TCanvas* ViewFEMesh::GetCanvas() { return m_canvas; }
73
74void ViewFEMesh::SetArea(double xmin, double ymin, double zmin, double xmax,
75 double ymax, double zmax) {
76
77 // Check range, assign if non-null
78 if (xmin == xmax || ymin == ymax) {
79 std::cout << m_className << "::SetArea:\n";
80 std::cout << " Null area range not permitted.\n";
81 return;
82 }
83 m_xMin = std::min(xmin, xmax);
84 m_yMin = std::min(ymin, ymax);
85 m_zMin = std::min(zmin, zmax);
86 m_xMax = std::max(xmin, xmax);
87 m_yMax = std::max(ymin, ymax);
88 m_zMax = std::max(zmin, zmax);
89
90 m_hasUserArea = true;
91}
92
93void ViewFEMesh::SetArea() { m_hasUserArea = false; }
94
95// The plotting functionality here is ported from Garfield
96// with some inclusion of code from ViewCell.cc
98
99 if (!m_component) {
100 std::cerr << m_className << "::Plot:\n";
101 std::cerr << " Component is not defined.\n";
102 return false;
103 }
104
105 double pmin = 0., pmax = 0.;
106 if (!m_component->GetVoltageRange(pmin, pmax)) {
107 std::cerr << m_className << "::Plot:\n";
108 std::cerr << " Component is not ready.\n";
109 return false;
110 }
111
112 // Get the bounding box.
113 if (!m_hasUserArea) {
114 std::cerr << m_className << "::Plot:\n";
115 std::cerr << " Bounding box cannot be determined.\n";
116 std::cerr << " Call SetArea first.\n";
117 return false;
118 }
119
120 // Set up a canvas if one does not already exist.
121 if (!m_canvas) {
122 m_canvas = new TCanvas();
123 m_canvas->SetTitle(m_label.c_str());
124 if (m_hasExternalCanvas) m_hasExternalCanvas = false;
125 }
126 m_canvas->Range(m_xMin, m_yMin, m_xMax, m_yMax);
127
128 // Plot the elements
129 ComponentCST* componentCST = dynamic_cast<ComponentCST*>(m_component);
130 if (componentCST) {
131 std::cout << m_className << "::Plot:\n";
132 std::cout << " The given component is a CST component.\n";
133 std::cout << " Method PlotCST is called now!.\n";
134 DrawCST(componentCST);
135 } else {
136 DrawElements();
137 }
138 m_canvas->Update();
139
140 return true;
141}
142
143// Set the projection plane: modified from ViewField.cc
144// to match functionality of Garfield
145void ViewFEMesh::SetPlane(const double fx, const double fy, const double fz,
146 const double x0, const double y0, const double z0) {
147
148 // Calculate 2 in-plane vectors for the normal vector
149 double fnorm = sqrt(fx * fx + fy * fy + fz * fz);
150 double dist = fx * x0 + fy * y0 + fz * z0;
151 if (fnorm > 0 && fx * fx + fz * fz > 0) {
152 project[0][0] = fz / sqrt(fx * fx + fz * fz);
153 project[0][1] = 0;
154 project[0][2] = -fx / sqrt(fx * fx + fz * fz);
155 project[1][0] = -fx * fy / (sqrt(fx * fx + fz * fz) * fnorm);
156 project[1][1] = (fx * fx + fz * fz) / (sqrt(fx * fx + fz * fz) * fnorm);
157 project[1][2] = -fy * fz / (sqrt(fx * fx + fz * fz) * fnorm);
158 project[2][0] = dist * fx / (fnorm * fnorm);
159 project[2][1] = dist * fy / (fnorm * fnorm);
160 project[2][2] = dist * fz / (fnorm * fnorm);
161 } else if (fnorm > 0 && fy * fy + fz * fz > 0) {
162 project[0][0] = (fy * fy + fz * fz) / (sqrt(fy * fy + fz * fz) * fnorm);
163 project[0][1] = -fx * fz / (sqrt(fy * fy + fz * fz) * fnorm);
164 project[0][2] = -fy * fz / (sqrt(fy * fy + fz * fz) * fnorm);
165 project[1][0] = 0;
166 project[1][1] = fz / sqrt(fy * fy + fz * fz);
167 project[1][2] = -fy / sqrt(fy * fy + fz * fz);
168 project[2][0] = dist * fx / (fnorm * fnorm);
169 project[2][1] = dist * fy / (fnorm * fnorm);
170 project[2][2] = dist * fz / (fnorm * fnorm);
171 } else {
172 std::cout << m_className << "::SetPlane:\n";
173 std::cout << " Normal vector has zero norm.\n";
174 std::cout << " No new projection set.\n";
175 }
176
177 // Store the plane description
178 plane[0] = fx;
179 plane[1] = fy;
180 plane[2] = fz;
181 plane[3] = dist;
182}
183
184// Set the default projection for the plane: copied from ViewField.cc
186
187 // Default projection: x-y at z=0
188 project[0][0] = 1;
189 project[1][0] = 0;
190 project[2][0] = 0;
191 project[0][1] = 0;
192 project[1][1] = 1;
193 project[2][1] = 0;
194 project[0][2] = 0;
195 project[1][2] = 0;
196 project[2][2] = 0;
197
198 // Plane description
199 plane[0] = 0;
200 plane[1] = 0;
201 plane[2] = 1;
202 plane[3] = 0;
203}
204
205// Set the x-axis.
206void ViewFEMesh::SetXaxis(TGaxis* ax) { m_xaxis = ax; }
207
208// Set the y-axis.
209void ViewFEMesh::SetYaxis(TGaxis* ay) { m_yaxis = ay; }
210
211// Set the x-axis title.
212void ViewFEMesh::SetXaxisTitle(const char* xtitle) {
213 m_axes->GetXaxis()->SetTitle(xtitle);
214}
215
216// Set the y-axis title.
217void ViewFEMesh::SetYaxisTitle(const char* ytitle) {
218 m_axes->GetYaxis()->SetTitle(ytitle);
219}
220
221// Create default axes
223
224 // Create a new x and y axis
225 m_xaxis = new TGaxis(m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
226 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
227 m_xMax - std::abs(m_xMax - m_xMin) * 0.1,
228 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
229 m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
230 m_xMax - std::abs(m_xMax - m_xMin) * 0.1, 2405, "x");
231 m_yaxis = new TGaxis(m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
232 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
233 m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
234 m_yMax - std::abs(m_yMax - m_yMin) * 0.1,
235 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
236 m_yMax - std::abs(m_yMax - m_yMin) * 0.1, 2405, "y");
237
238 // Label sizes
239 m_xaxis->SetLabelSize(0.025);
240 m_yaxis->SetLabelSize(0.025);
241
242 // Titles
243 m_xaxis->SetTitleSize(0.03);
244 m_xaxis->SetTitle("x [cm]");
245 m_yaxis->SetTitleSize(0.03);
246 m_yaxis->SetTitle("y [cm]");
247}
248
249// Use ROOT plotting functions to draw the mesh elements on the canvas.
250// General methodology ported from Garfield
251void ViewFEMesh::DrawElements() {
252
253 // Get the map boundaries from the component.
254 double mapxmax = m_component->mapxmax;
255 double mapxmin = m_component->mapxmin;
256 double mapymax = m_component->mapymax;
257 double mapymin = m_component->mapymin;
258 double mapzmax = m_component->mapzmax;
259 double mapzmin = m_component->mapzmin;
260
261 // Get the periodicities.
262 double sx = mapxmax - mapxmin;
263 double sy = mapymax - mapymin;
264 double sz = mapzmax - mapzmin;
265 const bool perX = m_component->m_xPeriodic || m_component->m_xMirrorPeriodic;
266 const bool perY = m_component->m_yPeriodic || m_component->m_yMirrorPeriodic;
267 const bool perZ = m_component->m_zPeriodic || m_component->m_zMirrorPeriodic;
268
269 // Clear the meshes and drift line lists.
270 m_mesh.clear();
271 m_driftLines.clear();
272
273 // Prepare the final projection matrix (the transpose of the 2D array
274 // "project").
275 double fnorm =
276 sqrt(plane[0] * plane[0] + plane[1] * plane[1] + plane[2] * plane[2]);
277 TArrayD dataProj(9);
278 dataProj[0] = project[0][0];
279 dataProj[1] = project[1][0];
280 dataProj[2] = plane[0] / fnorm;
281 dataProj[3] = project[0][1];
282 dataProj[4] = project[1][1];
283 dataProj[5] = plane[1] / fnorm;
284 dataProj[6] = project[0][2];
285 dataProj[7] = project[1][2];
286 dataProj[8] = plane[2] / fnorm;
287 TMatrixD projMat(3, 3, dataProj.GetArray());
288
289 // Calculate the determinant of the projection matrix.
290 double projDet =
291 projMat(0, 0) *
292 (projMat(1, 1) * projMat(2, 2) - projMat(1, 2) * projMat(2, 1)) -
293 projMat(0, 1) *
294 (projMat(1, 0) * projMat(2, 2) - projMat(1, 2) * projMat(2, 0)) +
295 projMat(0, 2) *
296 (projMat(1, 0) * projMat(2, 1) - projMat(1, 1) * projMat(2, 0));
297
298 // Calculate the inverse of the projection matrix for
299 // calculating coordinates in the viewing plane.
300 if (projDet != 0) {
301 projMat.Invert();
302 } else {
303 std::cerr << m_className << "::DrawElements:\n";
304 std::cerr << " Projection matrix is not invertible.\n";
305 std::cerr << " Finite element mesh will not be drawn.\n";
306 }
307
308 // Get the plane information.
309 double fx = plane[0];
310 double fy = plane[1];
311 double fz = plane[2];
312 double dist = plane[3];
313
314 // Construct two empty single-column matrices for use as coordinate vectors.
315 TMatrixD xMat(3, 1);
316
317 // Determine the number of periods present in the cell.
318 int nMaxX = 0, nMinX = 0;
319 int nMaxY = 0, nMinY = 0;
320 int nMaxZ = 0, nMinZ = 0;
321 if (perX) {
322 nMinX = int(m_xMin / sx) - 1;
323 nMaxX = int(m_xMax / sx) + 1;
324 }
325 if (perY) {
326 nMinY = int(m_yMin / sy) - 1;
327 nMaxY = int(m_yMax / sy) + 1;
328 }
329 if (perZ) {
330 nMinZ = int(m_zMin / sz) - 1;
331 nMaxZ = int(m_zMax / sz) + 1;
332 }
333
334 // Loop over all elements.
335 for (int elem = 0; elem < m_component->nElements; elem++) {
336
337 // Do not plot the drift medium.
338 if (m_component->materials[m_component->elements[elem].matmap]
339 .driftmedium &&
340 !(m_plotMeshBorders))
341 continue;
342 // Do not create Polygons for disabled materials
343 if (m_disabledMaterial[m_component->elements[elem].matmap]) continue;
344 // -- Tetrahedral elements
345
346 // Coordinates of vertices
347 double vx1, vy1, vz1;
348 double vx2, vy2, vz2;
349 double vx3, vy3, vz3;
350 double vx4, vy4, vz4;
351
352 // Get the color for this element (default to 1).
353 int colorID = m_colorMap.count(m_component->elements[elem].matmap);
354 if (colorID != 0)
355 colorID = m_colorMap[m_component->elements[elem].matmap];
356 else
357 colorID = 1;
358
359 // Get the fill color for this element (default colorID).
360 int colorID_fill =
361 m_colorMap_fill.count(m_component->elements[elem].matmap);
362 if (colorID_fill != 0)
363 colorID_fill = m_colorMap_fill[m_component->elements[elem].matmap];
364 else
365 colorID_fill = colorID;
366
367 // Loop over the periodicities in x.
368 for (int nx = nMinX; nx <= nMaxX; nx++) {
369
370 // Determine the x-coordinates of the tetrahedral vertices.
371 if (m_component->m_xMirrorPeriodic && nx != 2 * (nx / 2)) {
372 vx1 = mapxmin +
373 (mapxmax -
374 m_component->nodes[m_component->elements[elem].emap[0]].x) +
375 sx * nx;
376 vx2 = mapxmin +
377 (mapxmax -
378 m_component->nodes[m_component->elements[elem].emap[1]].x) +
379 sx * nx;
380 vx3 = mapxmin +
381 (mapxmax -
382 m_component->nodes[m_component->elements[elem].emap[2]].x) +
383 sx * nx;
384 vx4 = mapxmin +
385 (mapxmax -
386 m_component->nodes[m_component->elements[elem].emap[3]].x) +
387 sx * nx;
388 } else {
389 vx1 =
390 m_component->nodes[m_component->elements[elem].emap[0]].x + sx * nx;
391 vx2 =
392 m_component->nodes[m_component->elements[elem].emap[1]].x + sx * nx;
393 vx3 =
394 m_component->nodes[m_component->elements[elem].emap[2]].x + sx * nx;
395 vx4 =
396 m_component->nodes[m_component->elements[elem].emap[3]].x + sx * nx;
397 }
398
399 // Loop over the periodicities in y.
400 for (int ny = nMinY; ny <= nMaxY; ny++) {
401
402 // Determine the y-coordinates of the tetrahedral vertices.
403 if (m_component->m_yMirrorPeriodic && ny != 2 * (ny / 2)) {
404 vy1 = mapymin +
405 (mapymax -
406 m_component->nodes[m_component->elements[elem].emap[0]].y) +
407 sy * ny;
408 vy2 = mapymin +
409 (mapymax -
410 m_component->nodes[m_component->elements[elem].emap[1]].y) +
411 sy * ny;
412 vy3 = mapymin +
413 (mapymax -
414 m_component->nodes[m_component->elements[elem].emap[2]].y) +
415 sy * ny;
416 vy4 = mapymin +
417 (mapymax -
418 m_component->nodes[m_component->elements[elem].emap[3]].y) +
419 sy * ny;
420 } else {
421 vy1 = m_component->nodes[m_component->elements[elem].emap[0]].y +
422 sy * ny;
423 vy2 = m_component->nodes[m_component->elements[elem].emap[1]].y +
424 sy * ny;
425 vy3 = m_component->nodes[m_component->elements[elem].emap[2]].y +
426 sy * ny;
427 vy4 = m_component->nodes[m_component->elements[elem].emap[3]].y +
428 sy * ny;
429 }
430
431 // Loop over the periodicities in z.
432 for (int nz = nMinZ; nz <= nMaxZ; nz++) {
433
434 // Determine the z-coordinates of the tetrahedral vertices.
435 if (m_component->m_zMirrorPeriodic && nz != 2 * (nz / 2)) {
436 vz1 = mapzmin +
437 (mapzmax -
438 m_component->nodes[m_component->elements[elem].emap[0]].z) +
439 sz * nz;
440 vz2 = mapzmin +
441 (mapzmax -
442 m_component->nodes[m_component->elements[elem].emap[1]].z) +
443 sz * nz;
444 vz3 = mapzmin +
445 (mapzmax -
446 m_component->nodes[m_component->elements[elem].emap[2]].z) +
447 sz * nz;
448 vz4 = mapzmin +
449 (mapzmax -
450 m_component->nodes[m_component->elements[elem].emap[3]].z) +
451 sz * nz;
452 } else {
453 vz1 = m_component->nodes[m_component->elements[elem].emap[0]].z +
454 sz * nz;
455 vz2 = m_component->nodes[m_component->elements[elem].emap[1]].z +
456 sz * nz;
457 vz3 = m_component->nodes[m_component->elements[elem].emap[2]].z +
458 sz * nz;
459 vz4 = m_component->nodes[m_component->elements[elem].emap[3]].z +
460 sz * nz;
461 }
462
463 // Store the x and y coordinates of the relevant mesh vertices.
464 std::vector<double> vX;
465 std::vector<double> vY;
466
467 // Boolean values for determining what coordinates are in the viewing
468 // plane
469 bool in1 = false, in2 = false, in3 = false, in4 = false;
470
471 // Boolean value for determining whether a plane cut succeeded
472 bool planeCut = false;
473
474 // Value used to determine whether a vertex is in the plane.
475 double pcf = std::max(
476 std::max(std::max(std::max(std::max(std::max(std::abs(vx1),
477 std::abs(vy1)),
478 std::abs(vz1)),
479 std::abs(fx)),
480 std::abs(fy)),
481 std::abs(fz)),
482 std::abs(dist));
483
484 // First isolate the vertices that are in the viewing plane.
485 if (std::abs(fx * vx1 + fy * vy1 + fz * vz1 - dist) < 1.0e-4 * pcf) {
486 in1 = true;
487 }
488 if (std::abs(fx * vx2 + fy * vy2 + fz * vz2 - dist) < 1.0e-4 * pcf) {
489 in2 = true;
490 }
491 if (std::abs(fx * vx3 + fy * vy3 + fz * vz3 - dist) < 1.0e-4 * pcf) {
492 in3 = true;
493 }
494 if (std::abs(fx * vx4 + fy * vy4 + fz * vz4 - dist) < 1.0e-4 * pcf) {
495 in4 = true;
496 }
497
498 // Calculate the planar coordinates for those edges that are in the
499 // plane.
500 if (in1) {
501 PlaneCoords(vx1, vy1, vz1, projMat, xMat);
502 vX.push_back(xMat(0, 0));
503 vY.push_back(xMat(1, 0));
504 }
505 if (in2) {
506 PlaneCoords(vx2, vy2, vz2, projMat, xMat);
507 vX.push_back(xMat(0, 0));
508 vY.push_back(xMat(1, 0));
509 }
510 if (in3) {
511 PlaneCoords(vx3, vy3, vz3, projMat, xMat);
512 vX.push_back(xMat(0, 0));
513 vY.push_back(xMat(1, 0));
514 }
515 if (in4) {
516 PlaneCoords(vx4, vy4, vz4, projMat, xMat);
517 vX.push_back(xMat(0, 0));
518 vY.push_back(xMat(1, 0));
519 }
520
521 // Cut the sides that are not in the plane.
522 if (!(in1 || in2)) {
523 planeCut = PlaneCut(vx1, vy1, vz1, vx2, vy2, vz2, xMat);
524 if (planeCut) {
525 vX.push_back(xMat(0, 0));
526 vY.push_back(xMat(1, 0));
527 }
528 }
529 if (!(in1 || in3)) {
530 planeCut = PlaneCut(vx1, vy1, vz1, vx3, vy3, vz3, xMat);
531 if (planeCut) {
532 vX.push_back(xMat(0, 0));
533 vY.push_back(xMat(1, 0));
534 }
535 }
536 if (!(in1 || in4)) {
537 planeCut = PlaneCut(vx1, vy1, vz1, vx4, vy4, vz4, xMat);
538 if (planeCut) {
539 vX.push_back(xMat(0, 0));
540 vY.push_back(xMat(1, 0));
541 }
542 }
543 if (!(in2 || in3)) {
544 planeCut = PlaneCut(vx2, vy2, vz2, vx3, vy3, vz3, xMat);
545 if (planeCut) {
546 vX.push_back(xMat(0, 0));
547 vY.push_back(xMat(1, 0));
548 }
549 }
550 if (!(in2 || in4)) {
551 planeCut = PlaneCut(vx2, vy2, vz2, vx4, vy4, vz4, xMat);
552 if (planeCut) {
553 vX.push_back(xMat(0, 0));
554 vY.push_back(xMat(1, 0));
555 }
556 }
557 if (!(in3 || in4)) {
558 planeCut = PlaneCut(vx3, vy3, vz3, vx4, vy4, vz4, xMat);
559 if (planeCut) {
560 vX.push_back(xMat(0, 0));
561 vY.push_back(xMat(1, 0));
562 }
563 }
564
565 // Create a convex TPolyLine object connecting the points.
566 if (vX.size() >= 3) {
567
568 // Eliminate crossings of the polygon lines
569 // (known as "butterflies" in Garfield).
570 RemoveCrossings(vX, vY);
571
572 // Create vectors to store the clipped polygon.
573 std::vector<double> cX;
574 std::vector<double> cY;
575
576 // Clip the polygon to the view area.
577 ClipToView(vX, vY, cX, cY);
578
579 // If we have more than 2 vertices, add the polygon.
580 if (cX.size() > 2) {
581
582 // Again eliminate crossings of the polygon lines.
583 RemoveCrossings(cX, cY);
584
585 // Create the TPolyLine.
586 int polyPts = 0;
587 TPolyLine* poly = new TPolyLine();
588 poly->SetLineColor(colorID);
589 poly->SetFillColor(colorID_fill);
590 poly->SetLineWidth(3);
591
592 // Add all of the points.
593 for (int pt = 0; pt < (int)cX.size(); pt++) {
594
595 // Add this point.
596 poly->SetPoint(polyPts, cX[pt], cY[pt]);
597 polyPts++;
598 }
599
600 // Add the polygon to the mesh.
601 m_mesh.push_back(poly);
602 }
603 } // end TPolyLine construction if statement
604 } // end z-periodicity loop
605 } // end y-periodicity loop
606 } // end x-periodicity loop
607 } // end loop over elements
608
609 // If we have an associated ViewDrift, plot projections of the drift lines.
610 if (m_viewDrift) {
611 const int nDlines = m_viewDrift->m_driftLines.size();
612 for (int dline = 0; dline < nDlines; dline++) {
613 // Get the number of points.
614 const unsigned int npts = m_viewDrift->m_driftLines[dline].vect.size();
615 // Create a TPolyLine that is a 2D projection of the original.
616 TPolyLine* poly = new TPolyLine();
617 poly->SetLineColor(m_viewDrift->m_driftLines[dline].n);
618 int polyPts = 0;
619 for (unsigned int pt = 0; pt < npts; pt++) {
620 // Get the coordinates of this point.
621 double ptx = m_viewDrift->m_driftLines[dline].vect[pt].x;
622 double pty = m_viewDrift->m_driftLines[dline].vect[pt].y;
623 double ptz = m_viewDrift->m_driftLines[dline].vect[pt].z;
624 // Project this point onto the plane.
625 PlaneCoords(ptx, pty, ptz, projMat, xMat);
626 // Add this point if it is within the view.
627 if (xMat(0, 0) >= m_xMin && xMat(0, 0) <= m_xMax &&
628 xMat(1, 0) >= m_yMin && xMat(1, 0) <= m_yMax) {
629 poly->SetPoint(polyPts, xMat(0, 0), xMat(1, 0));
630 polyPts++;
631 }
632 } // end loop over points
633
634 // Add the drift line to the list.
635 m_driftLines.push_back(poly);
636
637 } // end loop over drift lines
638 } // end if(m_viewDrift != 0)
639
640 // Call the ROOT draw methods to plot the elements.
641 m_canvas->cd();
642
643 // Draw default axes by using a blank 2D histogram.
644 if (!m_xaxis && !m_yaxis && m_drawAxes) {
645 m_axes->GetXaxis()->SetLimits(m_xMin, m_xMax);
646 m_axes->GetYaxis()->SetLimits(m_yMin, m_yMax);
647 m_axes->Draw();
648 }
649
650 // Draw custom axes.
651 if (m_xaxis && m_drawAxes) m_xaxis->Draw();
652 if (m_yaxis && m_drawAxes) m_yaxis->Draw();
653
654 // Draw the mesh on the canvas.
655 for (int m = m_mesh.size(); m--;) {
656 if (m_plotMeshBorders || !m_fillMesh) m_mesh[m]->Draw("same");
657 if (m_fillMesh) m_mesh[m]->Draw("f:same");
658 }
659
660 // Draw the drift lines on the view.
661 for (int m = m_driftLines.size(); m--;) {
662 m_driftLines[m]->Draw("same");
663 }
664}
665
666void ViewFEMesh::DrawCST(ComponentCST* componentCST) {
667 /*The method is based on ViewFEMesh::Draw, thus the first part is copied from
668 * there.
669 * At the moment only x-y, x-z, and y-z are available due to the simple
670 * implementation.
671 * The advantage of this method is that there is no element loop and thus it
672 * is much
673 * faster.
674 */
675 // Get the map boundaries from the component
676 double mapxmax = m_component->mapxmax;
677 double mapxmin = m_component->mapxmin;
678 double mapymax = m_component->mapymax;
679 double mapymin = m_component->mapymin;
680 double mapzmax = m_component->mapzmax;
681 double mapzmin = m_component->mapzmin;
682
683 // Get the periodicities.
684 double sx = mapxmax - mapxmin;
685 double sy = mapymax - mapymin;
686 double sz = mapzmax - mapzmin;
687 const bool perX = m_component->m_xPeriodic || m_component->m_xMirrorPeriodic;
688 const bool perY = m_component->m_yPeriodic || m_component->m_yMirrorPeriodic;
689 const bool perZ = m_component->m_zPeriodic || m_component->m_zMirrorPeriodic;
690
691 // Clear the meshes and drift line lists
692 m_mesh.clear();
693 m_driftLines.clear();
694
695 // Prepare the final projection matrix (the transpose of the 2D array
696 // "project")
697 double fnorm =
698 sqrt(plane[0] * plane[0] + plane[1] * plane[1] + plane[2] * plane[2]);
699 TArrayD dataProj(9);
700 dataProj[0] = project[0][0];
701 dataProj[1] = project[1][0];
702 dataProj[2] = plane[0] / fnorm;
703 dataProj[3] = project[0][1];
704 dataProj[4] = project[1][1];
705 dataProj[5] = plane[1] / fnorm;
706 dataProj[6] = project[0][2];
707 dataProj[7] = project[1][2];
708 dataProj[8] = plane[2] / fnorm;
709 TMatrixD projMat(3, 3, dataProj.GetArray());
710
711 // Calculate the determinant of the projection matrix
712 double projDet =
713 projMat(0, 0) *
714 (projMat(1, 1) * projMat(2, 2) - projMat(1, 2) * projMat(2, 1)) -
715 projMat(0, 1) *
716 (projMat(1, 0) * projMat(2, 2) - projMat(1, 2) * projMat(2, 0)) +
717 projMat(0, 2) *
718 (projMat(1, 0) * projMat(2, 1) - projMat(1, 1) * projMat(2, 0));
719
720 // Calculate the inverse of the projection matrix for
721 // calculating coordinates in the viewing plane
722 if (projDet != 0) {
723 projMat.Invert();
724 } else {
725 std::cerr << m_className << "::DrawCST:\n";
726 std::cerr << " Projection matrix is not invertible.\n";
727 std::cerr << " Finite element mesh will not be drawn.\n";
728 }
729
730 // Construct two empty single-column matrices for use as coordinate vectors
731 TMatrixD xMat(3, 1);
732
733 // Determine the number of periods present in the cell.
734 int nMaxX = 0, nMinX = 0;
735 int nMaxY = 0, nMinY = 0;
736 int nMaxZ = 0, nMinZ = 0;
737 if (perX) {
738 nMinX = int(m_xMin / sx) - 1;
739 nMaxX = int(m_xMax / sx) + 1;
740 }
741 if (perY) {
742 nMinY = int(m_yMin / sy) - 1;
743 nMaxY = int(m_yMax / sy) + 1;
744 }
745 if (perZ) {
746 nMinZ = int(m_zMin / sz) - 1;
747 nMaxZ = int(m_zMax / sz) + 1;
748 }
749
750 int elem = 0;
751 std::vector<PolygonInfo> elements;
752 int nMinU = 0, nMaxU = 0, nMinV = 0, nMaxV = 0;
753 double mapumin = 0., mapumax = 0., mapvmin = 0., mapvmax = 0.;
754 double su = 0., sv = 0.;
755 bool mirroru = false, mirrorv = false;
756 double uMin, vMin, uMax, vMax;
757 unsigned int n_x, n_y, n_z;
758 componentCST->GetNumberOfMeshLines(n_x, n_y, n_z);
759 double e_xmin, e_xmax, e_ymin, e_ymax, e_zmin, e_zmax;
760 // xy view
761 if (plane[0] == 0 && plane[1] == 0 && plane[2] == 1) {
762 std::cout << m_className << "::DrawCST:\n";
763 std::cout << " Creating x-y mesh view.\n";
766 // calculate the z position
767 unsigned int i, j, z;
768 if (!componentCST->Coordinate2Index(0, 0, project[2][2], i, j, z)) {
769 std::cerr << "Could determine the position of the plane in z direction."
770 << std::endl;
771 return;
772 }
773 std::cout << " The plane position in z direction is: " << project[2][2]
774 << "\n";
775 nMinU = nMinX;
776 nMaxU = nMaxX;
777 nMinV = nMinY;
778 nMaxV = nMaxY;
779 uMin = m_xMin;
780 uMax = m_xMax;
781 vMin = m_yMin;
782 vMax = m_yMax;
783
784 mapumin = mapxmin;
785 mapumax = mapxmax;
786 mapvmin = mapymin;
787 mapvmax = mapymax;
788 su = sx;
789 sv = sy;
790 mirroru = perX;
791 mirrorv = perY;
792 for (unsigned int y = 0; y < (n_y - 1); y++) {
793 for (unsigned int x = 0; x < (n_x - 1); x++) {
794 elem = componentCST->Index2Element(x, y, z);
795 componentCST->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
796 e_zmin, e_zmax);
797 PolygonInfo tmp_info;
798 tmp_info.element = elem;
799 tmp_info.p1[0] = e_xmin;
800 tmp_info.p2[0] = e_xmax;
801 tmp_info.p3[0] = e_xmax;
802 tmp_info.p4[0] = e_xmin;
803 tmp_info.p1[1] = e_ymin;
804 tmp_info.p2[1] = e_ymin;
805 tmp_info.p3[1] = e_ymax;
806 tmp_info.p4[1] = e_ymax;
807 tmp_info.material = componentCST->GetElementMaterial(elem);
808 elements.push_back(tmp_info);
809 }
810 }
811 // xz-view
812 } else if (plane[0] == 0 && plane[1] == -1 && plane[2] == 0) {
813 std::cout << m_className << "::DrawCST:\n";
814 std::cout << " Creating x-z mesh view.\n";
817 // calculate the y position
818 unsigned int i = 0, j = 0, y = 0;
819 if (!componentCST->Coordinate2Index(0, project[2][1], 0, i, y, j)) {
820 std::cerr << "Could determine the position of the plane in y direction."
821 << std::endl;
822 return;
823 }
824 std::cout << " The plane position in y direction is: " << project[2][1]
825 << "\n";
826
827 nMinU = nMinX;
828 nMaxU = nMaxX;
829 nMinV = nMinZ;
830 nMaxV = nMaxZ;
831 uMin = m_xMin;
832 uMax = m_xMax;
833 vMin = m_zMin;
834 vMax = m_zMax;
835
836 mapumin = mapxmin;
837 mapumax = mapxmax;
838 mapvmin = mapzmin;
839 mapvmax = mapzmax;
840 su = sx;
841 sv = sz;
842 mirroru = perX;
843 mirrorv = perZ;
844 for (unsigned int z = 0; z < (n_z - 1); z++) {
845 for (unsigned int x = 0; x < (n_x - 1); x++) {
846 elem = componentCST->Index2Element(x, y, z);
847 componentCST->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
848 e_zmin, e_zmax);
849 PolygonInfo tmp_info;
850 tmp_info.element = elem;
851 tmp_info.p1[0] = e_xmin;
852 tmp_info.p2[0] = e_xmax;
853 tmp_info.p3[0] = e_xmax;
854 tmp_info.p4[0] = e_xmin;
855 tmp_info.p1[1] = e_zmin;
856 tmp_info.p2[1] = e_zmin;
857 tmp_info.p3[1] = e_zmax;
858 tmp_info.p4[1] = e_zmax;
859 tmp_info.material = componentCST->GetElementMaterial(elem);
860 elements.push_back(tmp_info);
861 }
862 }
863
864 // yz-view
865 } else if (plane[0] == -1 && plane[1] == 0 && plane[2] == 0) {
866 std::cout << m_className << "::DrawCST:\n";
867 std::cout << " Creating z-y mesh view.\n";
870 // calculate the x position
871 unsigned int i, j, x;
872 if (!componentCST->Coordinate2Index(project[2][0], 0, 0, x, i, j)) {
873 std::cerr << "Could determine the position of the plane in x direction."
874 << std::endl;
875 return;
876 }
877 std::cout << " The plane position in x direction is: " << project[2][0]
878 << "\n";
879 nMinU = nMinZ;
880 nMaxU = nMaxZ;
881 nMinV = nMinY;
882 nMaxV = nMaxY;
883 uMin = m_yMin;
884 uMax = m_yMax;
885 vMin = m_zMin;
886 vMax = m_zMax;
887
888 mapumin = mapzmin;
889 mapumax = mapzmax;
890 mapvmin = mapymin;
891 mapvmax = mapymax;
892 su = sz;
893 sv = sy;
894 mirroru = perZ;
895 mirrorv = perY;
896 for (unsigned int z = 0; z < (n_z - 1); z++) {
897 for (unsigned int y = 0; y < (n_y - 1); y++) {
898 elem = componentCST->Index2Element(x, y, z);
899 componentCST->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
900 e_zmin, e_zmax);
901 PolygonInfo tmp_info;
902 tmp_info.element = elem;
903 tmp_info.p1[0] = e_zmin;
904 tmp_info.p2[0] = e_zmax;
905 tmp_info.p3[0] = e_zmax;
906 tmp_info.p4[0] = e_zmin;
907 tmp_info.p1[1] = e_ymin;
908 tmp_info.p2[1] = e_ymin;
909 tmp_info.p3[1] = e_ymax;
910 tmp_info.p4[1] = e_ymax;
911 tmp_info.material = componentCST->GetElementMaterial(elem);
912 // Add the polygon to the mesh
913 elements.push_back(tmp_info);
914 }
915 }
916 } else {
917 std::cerr << m_className << "::DrawCST:\n";
918 std::cerr << " The given plane name is not known.\n";
919 std::cerr << " Please choose one of the following: xy, xz, yz.\n";
920 return;
921 }
922 std::cout << m_className << "::PlotCST:\n";
923 std::cout << " Number of elements in the projection of the unit cell:"
924 << elements.size() << std::endl;
925 std::vector<PolygonInfo>::iterator it;
926 std::vector<PolygonInfo>::iterator itend = elements.end();
927
928 for (int nu = nMinU; nu <= nMaxU; nu++) {
929 for (int nv = nMinV; nv <= nMaxV; nv++) {
930 it = elements.begin();
931 while (it != itend) {
932 if (m_disabledMaterial[(*it).material]) {
933 // do not create Polygons for disabled materials
934 it++;
935 continue;
936 }
937 int colorID = m_colorMap.count((*it).material);
938 if (colorID != 0)
939 colorID = m_colorMap[(*it).material];
940 else
941 colorID = 1;
942
943 // Get the fill color for this element (default colorID)
944 int colorID_fill = m_colorMap_fill.count((*it).material);
945 if (colorID_fill != 0)
946 colorID_fill = m_colorMap_fill[(*it).material];
947 else
948 colorID_fill = colorID;
949
950 TPolyLine* poly = new TPolyLine();
951 poly->SetLineColor(colorID);
952 poly->SetFillColor(colorID_fill);
953 if (m_plotMeshBorders)
954 poly->SetLineWidth(3);
955 else
956 poly->SetLineWidth(1);
957 // Add 4 points of the square
958 Double_t tmp_u[4], tmp_v[4];
959 if (mirroru && nu != 2 * (nu / 2)) {
960 // nu is odd
961 tmp_u[0] = mapumin + (mapumax - (*it).p1[0]) + su * nu;
962 tmp_u[1] = mapumin + (mapumax - (*it).p2[0]) + su * nu;
963 tmp_u[2] = mapumin + (mapumax - (*it).p3[0]) + su * nu;
964 tmp_u[3] = mapumin + (mapumax - (*it).p4[0]) + su * nu;
965 } else {
966 // nu is even
967 tmp_u[0] = (*it).p1[0] + su * nu;
968 tmp_u[1] = (*it).p2[0] + su * nu;
969 tmp_u[2] = (*it).p3[0] + su * nu;
970 tmp_u[3] = (*it).p4[0] + su * nu;
971 }
972 if (mirrorv && nv != 2 * (nv / 2)) {
973 tmp_v[0] = mapvmin + (mapvmax - (*it).p1[1]) + sv * nv;
974 tmp_v[1] = mapvmin + (mapvmax - (*it).p2[1]) + sv * nv;
975 tmp_v[2] = mapvmin + (mapvmax - (*it).p3[1]) + sv * nv;
976 tmp_v[3] = mapvmin + (mapvmax - (*it).p4[1]) + sv * nv;
977 } else {
978 tmp_v[0] = (*it).p1[1] + sv * nv;
979 tmp_v[1] = (*it).p2[1] + sv * nv;
980 tmp_v[2] = (*it).p3[1] + sv * nv;
981 tmp_v[3] = (*it).p4[1] + sv * nv;
982 }
983 if (tmp_u[0] < uMin || tmp_u[1] > uMax || tmp_v[0] < vMin ||
984 tmp_v[2] > vMax) {
985 it++;
986 continue;
987 }
988 poly->SetPoint(0, tmp_u[0], tmp_v[0]);
989 poly->SetPoint(1, tmp_u[1], tmp_v[1]);
990 poly->SetPoint(2, tmp_u[2], tmp_v[2]);
991 poly->SetPoint(3, tmp_u[3], tmp_v[3]);
992 // Add the polygon to the mesh
993 m_mesh.push_back(poly);
994 it++;
995 } // end element loop
996 } // end v-periodicity loop
997 } // end u-periodicity loop
998 std::cout << m_className << "::PlotCST:\n";
999 std::cout << " Number of polygons to be drawn:" << m_mesh.size()
1000 << std::endl;
1001 // If we have an associated ViewDrift, plot projections of the drift lines
1002 if (m_viewDrift) {
1003 const int nDlines = m_viewDrift->m_driftLines.size();
1004 for (int dline = 0; dline < nDlines; dline++) {
1005 // Get the number of points.
1006 const unsigned int npts = m_viewDrift->m_driftLines[dline].vect.size();
1007 // Create a TPolyLine that is a 2D projection of the original
1008 TPolyLine* poly = new TPolyLine();
1009 poly->SetLineColor(m_viewDrift->m_driftLines[dline].n);
1010 int polyPts = 0;
1011 for (unsigned int pt = 0; pt < npts; pt++) {
1012 // Get the coordinates of this point in the TPolyLine3D
1013 double ptx = m_viewDrift->m_driftLines[dline].vect[pt].x;
1014 double pty = m_viewDrift->m_driftLines[dline].vect[pt].y;
1015 double ptz = m_viewDrift->m_driftLines[dline].vect[pt].z;
1016 // Project this point onto the plane
1017 PlaneCoords(ptx, pty, ptz, projMat, xMat);
1018 // Add this point if it is within the view
1019 if (xMat(0, 0) >= uMin && xMat(0, 0) <= uMax && xMat(1, 0) >= vMin &&
1020 xMat(1, 0) <= vMax) {
1021 poly->SetPoint(polyPts, xMat(0, 0), xMat(1, 0));
1022 polyPts++;
1023 }
1024 } // end loop over points
1025 // Add the drift line to the list
1026 m_driftLines.push_back(poly);
1027 } // end loop over drift lines
1028 } // end if(m_viewDrift != 0)
1029
1030 // Call the ROOT draw methods to plot the elements
1031 m_canvas->cd();
1032 // Draw default axes by using a blank 2D histogram.
1033 if (!m_xaxis && !m_yaxis && m_drawAxes) {
1034 m_axes->GetXaxis()->SetLimits(uMin, uMax);
1035 m_axes->GetYaxis()->SetLimits(vMin, vMax);
1036 m_axes->Draw();
1037 }
1038 // Draw custom axes.
1039 if (m_xaxis && m_drawAxes) m_xaxis->Draw("");
1040 if (m_yaxis && m_drawAxes) m_yaxis->Draw("");
1041 // Draw the mesh on the canvas
1042 for (int m = m_mesh.size(); m--;) {
1043 if (m_plotMeshBorders || !m_fillMesh) m_mesh[m]->Draw("same");
1044 if (m_fillMesh) m_mesh[m]->Draw("f:sames");
1045 }
1046
1047 // Draw the drift lines on the view
1048 for (int m = m_driftLines.size(); m--;) {
1049 m_driftLines[m]->Draw("sames");
1050 }
1051 // TODO: Draw axes also at the end so that they are on top!
1052}
1053
1054// Returns true if the specified point is in the view region
1055bool ViewFEMesh::InView(double x, double y) {
1056 return (x >= m_xMin && x <= m_xMax && y >= m_yMin && y <= m_yMax);
1057}
1058
1059// Removes duplicate points and line crossings by correctly ordering
1060// the points in the provided vectors.
1061//
1062// NOTE: This is a 2D version of the BUTFLD method in Garfield. It
1063// follows the same general algorithm.
1064//
1065void ViewFEMesh::RemoveCrossings(std::vector<double>& x,
1066 std::vector<double>& y) {
1067
1068 // Determine element dimensions
1069 double xmin = x[0], xmax = x[0];
1070 double ymin = y[0], ymax = y[0];
1071 for (int i = 1; i < (int)x.size(); i++) {
1072 if (x[i] < xmin) xmin = x[i];
1073 if (x[i] > xmax) xmax = x[i];
1074 if (y[i] < ymin) ymin = y[i];
1075 if (y[i] > ymax) ymax = y[i];
1076 }
1077
1078 // First remove duplicate points
1079 double xtol = 1e-10 * std::abs(xmax - xmin);
1080 double ytol = 1e-10 * std::abs(ymax - ymin);
1081 for (int i = 0; i < (int)x.size(); i++) {
1082 for (int j = i + 1; j < (int)x.size(); j++) {
1083 if (std::abs(x[i] - x[j]) < xtol && std::abs(y[i] - y[j]) < ytol) {
1084 x.erase(x.begin() + j);
1085 y.erase(y.begin() + j);
1086 j--;
1087 }
1088 }
1089 }
1090
1091 // No crossings with 3 points or less
1092 if (x.size() <= 3) return;
1093
1094 // Save the polygon size so it is easily accessible
1095 int NN = x.size();
1096
1097 // Keep track of the number of attempts
1098 int attempts = 0;
1099
1100 // Exchange points until crossings are eliminated or we have attempted NN
1101 // times
1102 bool crossings = true;
1103 while (crossings && (attempts < NN)) {
1104
1105 // Assume we are done after this attempt.
1106 crossings = false;
1107
1108 for (int i = 1; i <= NN; i++) {
1109 for (int j = i + 2; j <= NN; j++) {
1110
1111 // End the j-loop if we have surpassed N and wrapped around to i
1112 if ((j + 1) > NN && 1 + (j % NN) >= i) break;
1113
1114 // Otherwise, detect crossings and attempt to eliminate them.
1115 else {
1116
1117 // Determine if we have a crossing.
1118 double xc = 0., yc = 0.;
1119 if (LinesCrossed(x[(i - 1) % NN], y[(i - 1) % NN], x[i % NN],
1120 y[i % NN], x[(j - 1) % NN], y[(j - 1) % NN],
1121 x[j % NN], y[j % NN], xc, yc)) {
1122
1123 // Swap each point from i towards j with each corresponding point
1124 // from j towards i.
1125 for (int k = 1; k <= (j - i) / 2; k++) {
1126
1127 double xs = x[(i + k - 1) % NN];
1128 double ys = y[(i + k - 1) % NN];
1129 x[(i + k - 1) % NN] = x[(j - k) % NN];
1130 y[(i + k - 1) % NN] = y[(j - k) % NN];
1131 x[(j - k) % NN] = xs;
1132 y[(j - k) % NN] = ys;
1133
1134 // Force another attempt
1135 crossings = true;
1136 }
1137 }
1138 }
1139 } // end loop over j
1140 } // end loop over i
1141
1142 // Increment the number of attempts
1143 attempts++;
1144
1145 } // end while(crossings)
1146
1147 if (attempts > NN) {
1148 std::cerr << m_className << "::RemoveCrossings:\n";
1149 std::cerr
1150 << " WARNING: Maximum attempts reached - crossings not removed.\n";
1151 }
1152}
1153
1154//
1155// Determines whether the line connecting points (x1,y1) and (x2,y2)
1156// and the line connecting points (u1,v1) and (u2,v2) cross somewhere
1157// between the 4 points. Sets the crossing point in (xc, yc).
1158//
1159// Ported from Garfield function CROSSD
1160//
1161bool ViewFEMesh::LinesCrossed(double x1, double y1, double x2, double y2,
1162 double u1, double v1, double u2, double v2,
1163 double& xc, double& yc) {
1164
1165 // Set the tolerances.
1166 double xtol =
1167 1.0e-10 *
1168 std::max(std::abs(x1),
1169 std::max(std::abs(x2), std::max(std::abs(u1), std::abs(u2))));
1170 double ytol =
1171 1.0e-10 *
1172 std::max(std::abs(y1),
1173 std::max(std::abs(y2), std::max(std::abs(v1), std::abs(v2))));
1174 if (xtol <= 0) xtol = 1.0e-10;
1175 if (ytol <= 0) ytol = 1.0e-10;
1176
1177 // Compute the distances and determinant (dx,dy) x (du,dv).
1178 double dy = y2 - y1;
1179 double dv = v2 - v1;
1180 double dx = x1 - x2;
1181 double du = u1 - u2;
1182 double det = dy * du - dx * dv;
1183
1184 // Check for crossing because one of the endpoints is on both lines.
1185 if (OnLine(x1, y1, x2, y2, u1, v1)) {
1186 xc = u1;
1187 yc = v1;
1188 return true;
1189 } else if (OnLine(x1, y1, x2, y2, u2, v2)) {
1190 xc = u2;
1191 yc = v2;
1192 return true;
1193 } else if (OnLine(u1, v1, u2, v2, x1, y1)) {
1194 xc = x1;
1195 yc = y1;
1196 return true;
1197 } else if (OnLine(u1, v1, u2, v2, x2, y2)) {
1198 xc = x2;
1199 yc = y2;
1200 return true;
1201 }
1202 // Check if the lines are parallel (zero determinant).
1203 else if (std::abs(det) < xtol * ytol)
1204 return false;
1205 // No special case: compute point of intersection.
1206 else {
1207
1208 // Solve crossing equations.
1209 xc = (du * (x1 * y2 - x2 * y1) - dx * (u1 * v2 - u2 * v1)) / det;
1210 yc = ((-1 * dv) * (x1 * y2 - x2 * y1) + dy * (u1 * v2 - u2 * v1)) / det;
1211
1212 // Determine if this point is on both lines.
1213 if (OnLine(x1, y1, x2, y2, xc, yc) && OnLine(u1, v1, u2, v2, xc, yc))
1214 return true;
1215 }
1216
1217 // The lines do not cross if we have reached this point.
1218 return false;
1219}
1220
1221//
1222// Determines whether the point (u,v) lies on the line connecting
1223// points (x1,y1) and (x2,y2).
1224//
1225// Ported from Garfield function ONLIND
1226//
1227bool ViewFEMesh::OnLine(double x1, double y1, double x2, double y2, double u,
1228 double v) {
1229
1230 // Set the tolerances
1231 double xtol =
1232 1.0e-10 * std::max(std::abs(x1), std::max(std::abs(x2), std::abs(u)));
1233 double ytol =
1234 1.0e-10 * std::max(std::abs(y1), std::max(std::abs(y2), std::abs(v)));
1235 if (xtol <= 0) xtol = 1.0e-10;
1236 if (ytol <= 0) ytol = 1.0e-10;
1237
1238 // To store the coordinates of the comparison point
1239 double xc = 0, yc = 0;
1240
1241 // Check if (u,v) is the point (x1,y1) or (x2,y2)
1242 if ((std::abs(x1 - u) <= xtol && std::abs(y1 - v) <= ytol) ||
1243 (std::abs(x2 - u) <= xtol && std::abs(y2 - v) <= ytol)) {
1244 return true;
1245 }
1246 // Check if the line is actually a point
1247 else if (std::abs(x1 - x2) <= xtol && std::abs(y1 - y2) <= ytol) {
1248 return false;
1249 }
1250 // Choose (x1,y1) as starting point if closer to (u,v)
1251 else if (std::abs(u - x1) + std::abs(v - y1) <
1252 std::abs(u - x2) + std::abs(v - y2)) {
1253
1254 // Compute the component of the line from (x1,y1) to (u,v)
1255 // along the line from (x1,y1) to (x2,y2)
1256 double dpar = ((u - x1) * (x2 - x1) + (v - y1) * (y2 - y1)) /
1257 ((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1));
1258
1259 // Determine the point on the line to which to compare (u,v)
1260 if (dpar < 0.0) {
1261 xc = x1;
1262 yc = y1;
1263 } else if (dpar > 1.0) {
1264 xc = x2;
1265 yc = y2;
1266 } else {
1267 xc = x1 + dpar * (x2 - x1);
1268 yc = y1 + dpar * (y2 - y1);
1269 }
1270 }
1271 // Choose (x2,y2) as starting point if closer to (u,v)
1272 else {
1273
1274 // Compute the component of the line from (x2,y2) to (u,v)
1275 // along the line from (x2,y2) to (x1,y1)
1276 double dpar = ((u - x2) * (x1 - x2) + (v - y2) * (y1 - y2)) /
1277 ((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1));
1278
1279 // Determine the point on the line to which to compare (u,v)
1280 if (dpar < 0.0) {
1281 xc = x2;
1282 yc = y2;
1283 } else if (dpar > 1.0) {
1284 xc = x1;
1285 yc = y1;
1286 } else {
1287 xc = x2 + dpar * (x1 - x2);
1288 yc = y2 + dpar * (y1 - y2);
1289 }
1290 }
1291
1292 // Compare the calculated point to (u,v)
1293 if (std::abs(u - xc) < xtol && std::abs(v - yc) < ytol) return true;
1294
1295 return false;
1296}
1297
1298// Ported from Garfield: determines the point of intersection, in planar
1299// coordinates, of a plane with the line connecting multiple points
1300// x1,y1,z1;x2,y2,z2: the world coordinates of the two points
1301// projMat;planeMat: the projection and plane matrices
1302// xMat: the resulting planar coordinates of the intersection point
1303bool ViewFEMesh::PlaneCut(double x1, double y1, double z1, double x2, double y2,
1304 double z2, TMatrixD& xMat) {
1305
1306 // Set up the matrix for cutting edges not in the plane
1307 TArrayD dataCut(9);
1308 TMatrixD cutMat(3, 3);
1309 dataCut[0] = project[0][0];
1310 dataCut[1] = project[1][0];
1311 dataCut[2] = x1 - x2;
1312 dataCut[3] = project[0][1];
1313 dataCut[4] = project[1][1];
1314 dataCut[5] = y1 - y2;
1315 dataCut[6] = project[0][2];
1316 dataCut[7] = project[1][2];
1317 dataCut[8] = z1 - z2;
1318 cutMat.SetMatrixArray(dataCut.GetArray());
1319
1320 // Calculate the determinant of the cut matrix
1321 double cutDet =
1322 cutMat(0, 0) *
1323 (cutMat(1, 1) * cutMat(2, 2) - cutMat(1, 2) * cutMat(2, 1)) -
1324 cutMat(0, 1) *
1325 (cutMat(1, 0) * cutMat(2, 2) - cutMat(1, 2) * cutMat(2, 0)) +
1326 cutMat(0, 2) *
1327 (cutMat(1, 0) * cutMat(2, 1) - cutMat(1, 1) * cutMat(2, 0));
1328
1329 // Do not proceed if the matrix is singular
1330 if (cutDet == 0) return false;
1331
1332 // Set up a coordinate vector (RHS of equation)
1333 TArrayD dataCoords(3);
1334 TMatrixD coordMat(3, 1);
1335 dataCoords[0] = x1 - project[2][0];
1336 dataCoords[1] = y1 - project[2][1];
1337 dataCoords[2] = z1 - project[2][2];
1338 coordMat.SetMatrixArray(dataCoords.GetArray());
1339
1340 // Invert the cut matrix and multiply to get the solution
1341 cutMat.Invert();
1342 xMat = cutMat * coordMat;
1343
1344 // Return success if the plane point is between the two vertices
1345 if (xMat(2, 0) < 0 || xMat(2, 0) > 1) return false;
1346 return true;
1347}
1348
1349// Ported from Garfield: calculates the planar coordinates
1350// x,y,z: original world coordinates
1351// projMat: the projection matrix
1352// xMat: the resulting planar coordinates in single-column (vector) form
1353bool ViewFEMesh::PlaneCoords(double x, double y, double z,
1354 const TMatrixD& projMat, TMatrixD& xMat) {
1355
1356 // Set up the coordinate vector
1357 TArrayD dataCoords(3);
1358 TMatrixD coordMat(3, 1);
1359 dataCoords[0] = x;
1360 dataCoords[1] = y;
1361 dataCoords[2] = z;
1362 coordMat.SetMatrixArray(dataCoords.GetArray());
1363 xMat = projMat * coordMat;
1364
1365 return true;
1366}
1367
1368// Ported from Garfield (function INTERD):
1369// Returns true if the point (x,y) is inside of the specified polygon.
1370// x: the x-coordinate
1371// y: the y-coordinate
1372// px: the x-vertices of the polygon
1373// py: the y-vertices of the polygon
1374// edge: a variable set to true if the point is located on the polygon edge
1375bool ViewFEMesh::IsInPolygon(double x, double y, std::vector<double>& px,
1376 std::vector<double>& py, bool& edge) {
1377
1378 // Get the number and coordinates of the polygon vertices.
1379 int pN = (int)px.size();
1380
1381 // Handle the special case of less than 2 vertices.
1382 if (pN < 2) return false;
1383 // Handle the special case of exactly 2 vertices (a line).
1384 if (pN == 2) return OnLine(px[0], py[0], px[1], py[1], x, y);
1385
1386 // Set the minimum and maximum coordinates of all polygon vertices.
1387 double px_min = px[0], py_min = py[0];
1388 double px_max = px[0], py_max = py[0];
1389 for (int i = 0; i < pN; i++) {
1390 px_min = std::min(px_min, px[i]);
1391 py_min = std::min(py_min, py[i]);
1392 px_max = std::max(px_max, px[i]);
1393 py_max = std::max(py_max, py[i]);
1394 }
1395
1396 // Set the tolerances
1397 double xtol = 1.0e-10 * std::max(std::abs(px_min), std::abs(px_max));
1398 double ytol = 1.0e-10 * std::max(std::abs(py_min), std::abs(py_max));
1399 if (xtol <= 0) xtol = 1.0e-10;
1400 if (ytol <= 0) ytol = 1.0e-10;
1401
1402 // If we have essentially one x value, check to see if y is in range.
1403 if (std::abs(px_max - px_min) < xtol) {
1404 edge = (y > (py_min - ytol) && y < (py_max + ytol) &&
1405 std::abs(px_max + px_min - 2 * x) < xtol);
1406 return false;
1407 }
1408 // If we have essentially one y value, check to see if x is in range.
1409 if (std::abs(py_max - py_min) < ytol) {
1410 edge = (x > (px_min - xtol) && x < (px_max + xtol) &&
1411 std::abs(py_max + py_min - 2 * y) < ytol);
1412 return false;
1413 }
1414
1415 // Set "infinity" points.
1416 double xinf = px_min - std::abs(px_max - px_min);
1417 double yinf = py_min - std::abs(py_max - py_min);
1418
1419 // Loop until successful or maximum iterations (100) reached.
1420 int niter = 0;
1421 bool done = false;
1422 int ncross = 0;
1423
1424 while (!done && niter < 100) {
1425
1426 // Assume we will finish on this loop.
1427 done = true;
1428
1429 // Loop over all edges, counting the number of edges crossed by a line
1430 // extending from (x,y) to (xinf,yinf).
1431 ncross = 0;
1432 for (int i = 0; (done && i < pN); i++) {
1433
1434 // Determine whether the point lies on the edge.
1435 if (OnLine(px[i % pN], py[i % pN], px[(i + 1) % pN], py[(i + 1) % pN], x,
1436 y)) {
1437
1438 edge = true;
1439 return false;
1440 }
1441
1442 // Determine whether this edge is crossed; if so increment the counter.
1443 double xc = 0., yc = 0.;
1444 if (LinesCrossed(x, y, xinf, yinf, px[i % pN], py[i % pN],
1445 px[(i + 1) % pN], py[(i + 1) % pN], xc, yc))
1446 ncross++;
1447
1448 // Ensure this vertex is not crossed by the line from (x,y)
1449 // to (xinf,yinf); if so recompute (xinf,yinf) and start over.
1450 if (OnLine(x, y, xinf, yinf, px[i], py[i])) {
1451
1452 // Recompute (xinf,yinf).
1453 xinf = px_min - RndmUniform() * std::abs(px_max - xinf);
1454 yinf = py_min - RndmUniform() * std::abs(py_max - yinf);
1455
1456 // Start over.
1457 done = false;
1458 niter++;
1459 }
1460 }
1461 }
1462
1463 // If we failed to finish iterating, return false.
1464 if (niter >= 100) {
1465 std::cerr << m_className << "::IsInPolygon: unable to determine whether ("
1466 << x << ", " << y << ") is inside a polygon. Returning false.\n";
1467 return false;
1468 }
1469
1470 // Point is inside for an odd, nonzero number of crossings.
1471 return (ncross != 2 * (ncross / 2));
1472}
1473
1474// Ported from Garfield (method GRCONV):
1475// Clip the specified polygon to the view region; return the clipped polygon.
1476// px: the x-vertices of the polygon
1477// py: the y-vertices of the polygon
1478// cx: to contain the x-vertices of the clipped polygon
1479// cy: to contain the y-vertices of the clipped polygon
1480void ViewFEMesh::ClipToView(std::vector<double>& px, std::vector<double>& py,
1481 std::vector<double>& cx, std::vector<double>& cy) {
1482
1483 // Get the number and coordinates of the polygon vertices.
1484 int pN = (int)px.size();
1485
1486 // Clear the vectors to contain the final polygon.
1487 cx.clear();
1488 cy.clear();
1489
1490 // Set up the view vertices (counter-clockwise, starting at upper left).
1491 std::vector<double> vx;
1492 vx.push_back(m_xMin);
1493 vx.push_back(m_xMax);
1494 vx.push_back(m_xMax);
1495 vx.push_back(m_xMin);
1496 std::vector<double> vy;
1497 vy.push_back(m_yMax);
1498 vy.push_back(m_yMax);
1499 vy.push_back(m_yMin);
1500 vy.push_back(m_yMin);
1501 int vN = (int)vx.size();
1502
1503 // Do nothing if we have less than 2 points.
1504 if (pN < 2) return;
1505
1506 // Loop over the polygon vertices.
1507 for (int i = 0; i < pN; i++) {
1508
1509 // Flag for skipping check for edge intersection.
1510 bool skip = false;
1511
1512 // Loop over the view vertices.
1513 for (int j = 0; j < vN; j++) {
1514
1515 // Determine whether this vertex lies on a view edge:
1516 // if so add the vertex to the final polygon.
1517 if (OnLine(vx[j % vN], vy[j % vN], vx[(j + 1) % vN], vy[(j + 1) % vN],
1518 px[i], py[i])) {
1519
1520 // Add the vertex.
1521 cx.push_back(px[i]);
1522 cy.push_back(py[i]);
1523
1524 // Skip edge intersection check in this case.
1525 skip = true;
1526 }
1527
1528 // Determine whether a corner of the view area lies on this edge:
1529 // if so add the corner to the final polygon.
1530 if (OnLine(px[i % pN], py[i % pN], px[(i + 1) % pN], py[(i + 1) % pN],
1531 vx[j], vy[j])) {
1532
1533 // Add the vertex.
1534 cx.push_back(vx[j]);
1535 cy.push_back(vy[j]);
1536
1537 // Skip edge intersection check in this case.
1538 skip = true;
1539 }
1540 }
1541
1542 // If we have not skipped the edge intersection check, look for an
1543 // intersection between this edge and the view edges.
1544 if (!skip) {
1545
1546 // Loop over the view vertices.
1547 for (int j = 0; j < vN; j++) {
1548
1549 // Check for a crossing with this edge;
1550 // if one exists, add the crossing point.
1551 double xc = 0., yc = 0.;
1552 if (LinesCrossed(vx[j % vN], vy[j % vN], vx[(j + 1) % vN],
1553 vy[(j + 1) % vN], px[i % pN], py[i % pN],
1554 px[(i + 1) % pN], py[(i + 1) % pN], xc, yc)) {
1555
1556 // Add a vertex.
1557 cx.push_back(xc);
1558 cy.push_back(yc);
1559 }
1560 }
1561 }
1562 }
1563
1564 // Find all view field vertices inside the polygon.
1565 for (int j = 0; j < vN; j++) {
1566
1567 // Test whether this vertex is inside the polygon.
1568 // If so, add it to the final polygon.
1569 bool edge = false;
1570 if (IsInPolygon(vx[j], vy[j], px, py, edge)) {
1571
1572 // Add the view vertex.
1573 cx.push_back(vx[j]);
1574 cy.push_back(vy[j]);
1575 }
1576 }
1577
1578 // Find all polygon vertices inside the box.
1579 for (int i = 0; i < pN; i++) {
1580
1581 // Test whether this vertex is inside the view.
1582 // If so, add it to the final polygon.
1583 bool edge = false;
1584 if (IsInPolygon(px[i], py[i], vx, vy, edge)) {
1585
1586 // Add the polygon vertex.
1587 cx.push_back(px[i]);
1588 cy.push_back(py[i]);
1589 }
1590 }
1591}
1592}
bool m_zMirrorPeriodic
Mirror periodicity in z.
bool m_yPeriodic
Simple periodicity in y.
bool m_yMirrorPeriodic
Mirror periodicity in y.
bool m_xPeriodic
Simple periodicity in x.
bool m_zPeriodic
Simple periodicity in z.
bool m_xMirrorPeriodic
Mirror periodicity in x.
Base class for components based on finite-element field maps.
std::vector< Material > materials
virtual bool GetVoltageRange(double &vmin, double &vmax)
Calculate the voltage range [V].
std::vector< Element > elements
void SetXaxisTitle(const char *xtitle)
Definition: ViewFEMesh.cc:212
void CreateDefaultAxes()
Create a default set of custom-made axes.
Definition: ViewFEMesh.cc:222
void SetCanvas(TCanvas *c)
Set the canvas to be painted on.
Definition: ViewFEMesh.cc:61
ViewFEMesh()
Constructor.
Definition: ViewFEMesh.cc:13
void SetComponent(ComponentFieldMap *comp)
Set the component from which to retrieve the mesh and field.
Definition: ViewFEMesh.cc:50
TCanvas * GetCanvas()
Definition: ViewFEMesh.cc:72
void SetYaxisTitle(const char *ytitle)
Definition: ViewFEMesh.cc:217
bool Plot()
Plot method to be called by user.
Definition: ViewFEMesh.cc:97
void SetPlane(double fx, double fy, double fz, double x0, double y0, double z0)
Set the projection plane.
Definition: ViewFEMesh.cc:145
void SetArea()
Set area to be plotted to the default.
Definition: ViewFEMesh.cc:93
~ViewFEMesh()
Destructor.
Definition: ViewFEMesh.cc:45
void SetXaxis(TGaxis *ax)
Definition: ViewFEMesh.cc:206
void SetYaxis(TGaxis *ay)
Definition: ViewFEMesh.cc:209
void SetDefaultProjection()
Reset the projection plane.
Definition: ViewFEMesh.cc:185
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
PlottingEngineRoot plottingEngine
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314