Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ComponentNeBem3d.cc
Go to the documentation of this file.
1#include <algorithm>
2#include <cmath>
3#include <fstream>
4#include <iomanip>
5#include <iostream>
6#include <numeric>
7#include <set>
8#include <vector>
9#include <cfloat>
10
11#include "neBEMInterface.h"
12#include "neBEM.h"
13#include "NR.h"
14
18#include "Garfield/Polygon.hh"
19
20namespace {
21
22unsigned int NextPoint(const unsigned int i, const unsigned int n) {
23 const unsigned int j = i + 1;
24 return j < n ? j : 0;
25}
26
27unsigned int PrevPoint(const unsigned int i, const unsigned int n) {
28 return i > 0 ? i - 1 : n - 1;
29}
30
31double Mag(const std::array<double, 3>& a) {
32 return sqrt(a[0] * a[0] + a[1] * a[1] + a[2] * a[2]);
33}
34
35std::array<double, 3> UnitVector(const std::array<double, 3>& a) {
36
37 const double mag = Mag(a);
38 if (mag < 1.e-12) return a;
39 const std::array<double, 3> b = {a[0] / mag, a[1] / mag, a[2] / mag};
40 return b;
41}
42std::array<double, 3> CrossProduct(const std::array<double, 3>& u,
43 const std::array<double, 3>& v) {
44
45 const std::array<double, 3> w = {u[1] * v[2] - u[2] * v[1],
46 u[2] * v[0] - u[0] * v[2],
47 u[0] * v[1] - u[1] * v[0]};
48 return w;
49}
50
51std::array<double, 3> LocalToGlobal(const double x, const double y,
52 const double z,
53 const std::array<std::array<double, 3>, 3>& dcos,
54 const std::array<double, 3>& t) {
55
56 // Initial vector.
57 std::array<double, 4> u = {x, y, z, 1.};
58
59 std::array<std::array<double, 4>, 4> rot;
60 rot[0] = {dcos[0][0], dcos[1][0], dcos[2][0], 0.};
61 rot[1] = {dcos[0][1], dcos[1][1], dcos[2][1], 0.};
62 rot[2] = {dcos[0][2], dcos[1][2], dcos[2][2], 0.};
63 rot[3] = {0., 0., 0., 1.};
64
65 std::array<double, 4> v = {0., 0., 0., 0.};
66 for (int i = 0; i < 4; ++i) {
67 for (int j = 0; j < 4; ++j) {
68 v[i] += rot[i][j] * u[j];
69 }
70 }
71 const std::array<double, 3> a = {v[0] + t[0], v[1] + t[1], v[2] + t[2]};
72 return a;
73}
74
75/// Compute lambda for a point on a line (0 = start, 1 = end).
76double Lambda(const double x1, const double x0, const double x2,
77 const double y1, const double y0, const double y2) {
78 // Segment of zero length.
79 if ((x1 - x2) == 0. && (y1 - y2) == 0.) {
80 std::cerr << "ComponentNeBem3d::Lambda: Zero length segment.\n";
81 return 2.;
82 }
83
84 double xl = 0.;
85 const double dx1 = x0 - x1;
86 const double dy1 = y0 - y1;
87 const double dx2 = x0 - x2;
88 const double dy2 = y0 - y2;
89 if (dx1 * dx1 + dy1 * dy1 < dx2 * dx2 + dy2 * dy2) {
90 // Point nearer to (x1, y1).
91 if (fabs(y1 - y2) > fabs(x1 - x2)) {
92 xl = dy1 / (y2 - y1);
93 } else {
94 xl = dx1 / (x2 - x1);
95 }
96 } else {
97 // Point nearer to (x2, y2).
98 if (fabs(y1 - y2) > fabs(x1 - x2)) {
99 xl = 1. - dy2 / (y1 - y2);
100 } else {
101 xl = 1. - dx2 / (x1 - x2);
102 }
103 }
104 return xl;
105}
106
107/// Determine whether a point (u, v) lies on a straight line
108/// (x1, y1) to (x2, y2).
109bool OnLine(const double x1, const double y1, const double x2, const double y2,
110 const double u, const double v) {
111 // Set tolerances.
112 double epsx = 1.e-10 * std::max({fabs(x1), fabs(x2), fabs(u)});
113 double epsy = 1.e-10 * std::max({fabs(y1), fabs(y2), fabs(v)});
114 epsx = std::max(1.e-10, epsx);
115 epsy = std::max(1.e-10, epsy);
116
117 if ((fabs(x1 - u) <= epsx && fabs(y1 - v) <= epsy) ||
118 (fabs(x2 - u) <= epsx && fabs(y2 - v) <= epsy)) {
119 // Point to be examined coincides with start or end.
120 return true;
121 } else if (fabs(x1 - x2) <= epsx && fabs(y1 - y2) <= epsy) {
122 // The line (x1, y1) to (x2, y2) is in fact a point.
123 return false;
124 }
125 double xc = 0., yc = 0.;
126 if (fabs(u - x1) + fabs(v - y1) < fabs(u - x2) + fabs(v - y2)) {
127 // (u, v) is nearer to (x1, y1).
128 const double dx = (x2 - x1);
129 const double dy = (y2 - y1);
130 const double xl = ((u - x1) * dx + (v - y1) * dy) / (dx * dx + dy * dy);
131 if (xl < 0.) {
132 xc = x1;
133 yc = y1;
134 } else if (xl > 1.) {
135 xc = x2;
136 yc = y2;
137 } else {
138 xc = x1 + xl * dx;
139 yc = y1 + xl * dy;
140 }
141 } else {
142 // (u, v) is nearer to (x2, y2).
143 const double dx = (x1 - x2);
144 const double dy = (y1 - y2);
145 const double xl = ((u - x2) * dx + (v - y2) * dy) / (dx * dx + dy * dy);
146 if (xl < 0.) {
147 xc = x2;
148 yc = y2;
149 } else if (xl > 1.) {
150 xc = x1;
151 yc = y1;
152 } else {
153 xc = x2 + xl * dx;
154 yc = y2 + xl * dy;
155 }
156 }
157 // See whether the point is on the line.
158 if (fabs(u - xc) < epsx && fabs(v - yc) < epsy) {
159 return true;
160 }
161 return false;
162}
163
164/// Determine whether the 2 straight lines (x1, y1) to (x2, y2)
165/// and (u1, v1) to (u2, v2) cross at an intermediate point for both lines.
166bool Crossing(const double x1, const double y1, const double x2,
167 const double y2, const double u1, const double v1,
168 const double u2, const double v2, double& xc, double& yc) {
169 /// Matrix to compute the crossing point.
170 std::array<std::array<double, 2>, 2> a;
171 a[0][0] = y2 - y1;
172 a[0][1] = v2 - v1;
173 a[1][0] = x1 - x2;
174 a[1][1] = u1 - u2;
175 const double det = a[0][0] * a[1][1] - a[1][0] * a[0][1];
176 // Initial values.
177 xc = 0.;
178 yc = 0.;
179 // Set tolerances.
180 double epsx = 1.e-10 * std::max({fabs(x1), fabs(x2), fabs(u1), fabs(u2)});
181 double epsy = 1.e-10 * std::max({fabs(y1), fabs(y2), fabs(v1), fabs(v2)});
182 epsx = std::max(epsx, 1.e-10);
183 epsy = std::max(epsy, 1.e-10);
184 // Check for a point of one line located on the other line.
185 if (OnLine(x1, y1, x2, y2, u1, v1)) {
186 xc = u1;
187 yc = v1;
188 return true;
189 } else if (OnLine(x1, y1, x2, y2, u2, v2)) {
190 xc = u2;
191 yc = v2;
192 return true;
193 } else if (OnLine(u1, v1, u2, v2, x1, y1)) {
194 xc = x1;
195 yc = y1;
196 return true;
197 } else if (OnLine(u1, v1, u2, v2, x2, y2)) {
198 xc = x2;
199 yc = y2;
200 return true;
201 } else if (fabs(det) < epsx * epsy) {
202 // Parallel, non-touching.
203 return false;
204 }
205 // Crossing, non-trivial lines: solve crossing equations.
206 const double aux = a[1][1];
207 a[1][1] = a[0][0] / det;
208 a[0][0] = aux / det;
209 a[1][0] = -a[1][0] / det;
210 a[0][1] = -a[0][1] / det;
211 // Compute crossing point.
212 xc = a[0][0] * (x1 * y2 - x2 * y1) + a[1][0] * (u1 * v2 - u2 * v1);
213 yc = a[0][1] * (x1 * y2 - x2 * y1) + a[1][1] * (u1 * v2 - u2 * v1);
214 // See whether the crossing point is on both lines.
215 if (OnLine(x1, y1, x2, y2, xc, yc) && OnLine(u1, v1, u2, v2, xc, yc)) {
216 // Intersecting lines.
217 return true;
218 }
219 // Crossing point not on both lines.
220 return false;
221}
222
223void AddPoints(const std::vector<double>& xp1, const std::vector<double>& yp1,
224 const std::vector<double>& xp2, const std::vector<double>& yp2,
225 std::vector<double>& xl, std::vector<double>& yl,
226 std::vector<int>& flags, std::vector<double>& qs,
227 const double epsx, const double epsy) {
228 struct Point {
229 double x;
230 double y;
231 int flag;
232 double q;
233 };
234
235 std::vector<Point> points;
236
237 const unsigned int np1 = xp1.size();
238 const unsigned int np2 = xp2.size();
239 for (unsigned int i = 0; i < np1; ++i) {
240 const double xi0 = xp1[i];
241 const double yi0 = yp1[i];
242 const double xi1 = xp1[NextPoint(i, np1)];
243 const double yi1 = yp1[NextPoint(i, np1)];
244 // Add the vertex.
245 Point p1;
246 p1.x = xi0;
247 p1.y = yi0;
248 p1.flag = 1;
249 p1.q = 0.;
250 // If also on 2 or vertex of 2, flag it as crossing or foreign.
251 for (unsigned int j = 0; j < np2; ++j) {
252 const double xj0 = xp2[j];
253 const double yj0 = yp2[j];
254 if (fabs(xj0 - xi0) < epsx && fabs(yj0 - yi0) < epsy) {
255 p1.flag = 2;
256 }
257 const double xj1 = xp2[NextPoint(j, np2)];
258 const double yj1 = yp2[NextPoint(j, np2)];
259 if (OnLine(xj0, yj0, xj1, yj1, xi0, yi0) &&
260 (fabs(xj0 - xi0) > epsx || fabs(yj0 - yi0) > epsy) &&
261 (fabs(xj1 - xi0) > epsx || fabs(yj1 - yi0) > epsy)) {
262 p1.flag = 3;
263 }
264 }
265 points.push_back(std::move(p1));
266 // Go over the line segments of the other polygon.
267 std::vector<Point> pointsOther;
268 for (unsigned int j = 0; j < np2; ++j) {
269 const double xj0 = xp2[j];
270 const double yj0 = yp2[j];
271 // Add vertices of 2 that are on this line.
272 if (OnLine(xi0, yi0, xi1, yi1, xj0, yj0) &&
273 (fabs(xi0 - xj0) > epsx || fabs(yi0 - yj0) > epsy) &&
274 (fabs(xi1 - xj0) > epsx || fabs(yi1 - yj0) > epsy)) {
275 Point p2;
276 p2.x = xj0;
277 p2.y = yj0;
278 p2.flag = 2;
279 pointsOther.push_back(std::move(p2));
280 }
281 const double xj1 = xp2[NextPoint(j, np2)];
282 const double yj1 = yp2[NextPoint(j, np2)];
283 // Add crossing points.
284 double xc = 0., yc = 0.;
285 bool add = Crossing(xi0, yi0, xi1, yi1, xj0, yj0, xj1, yj1, xc, yc);
286 if (add) {
287 if ((fabs(xi0 - xc) < epsx && fabs(yi0 - yc) < epsy) ||
288 (fabs(xi1 - xc) < epsx && fabs(yi1 - yc) < epsy) ||
289 (fabs(xj0 - xc) < epsx && fabs(yj0 - yc) < epsy) ||
290 (fabs(xj1 - xc) < epsx && fabs(yj1 - yc) < epsy)) {
291 add = false;
292 }
293 if ((fabs(xi0 - xj0) < epsx && fabs(yi0 - yj0) < epsy) ||
294 (fabs(xi0 - xj1) < epsx && fabs(yi0 - yj1) < epsy) ||
295 (fabs(xi1 - xj0) < epsx && fabs(yi1 - yj0) < epsy) ||
296 (fabs(xi1 - xj1) < epsx && fabs(yi1 - yj1) < epsy)) {
297 add = false;
298 }
299 }
300 if (add) {
301 Point p2;
302 p2.x = xc;
303 p2.y = yc;
304 p2.flag = 3;
305 pointsOther.push_back(std::move(p2));
306 }
307 }
308 // Compute the lambdas for these points.
309 for (auto& p : pointsOther) {
310 p.q = Lambda(xi0, p.x, xi1, yi0, p.y, yi1);
311 }
312 // Sort the list by using the lambdas.
313 std::sort(
314 pointsOther.begin(), pointsOther.end(),
315 [](const Point& lhs, const Point& rhs) { return (lhs.q < rhs.q); });
316 points.insert(points.end(), pointsOther.begin(), pointsOther.end());
317 }
318
319 for (const Point& p : points) {
320 xl.push_back(p.x);
321 yl.push_back(p.y);
322 flags.push_back(p.flag);
323 qs.push_back(p.q);
324 }
325}
326
327/// Determine whether 2 panels are equal.
328bool Equal(const Garfield::Panel& panel1, const Garfield::Panel& panel2,
329 const double epsx, const double epsy) {
330 const auto& xp1 = panel1.xv;
331 const auto& yp1 = panel1.yv;
332 const auto& xp2 = panel2.xv;
333 const auto& yp2 = panel2.yv;
334 if (xp1.empty() || xp2.empty()) return false;
335 const unsigned int np1 = xp1.size();
336 const unsigned int np2 = xp2.size();
337
338 // Compare all points of 1 with all points of 2.
339 for (unsigned int i = 0; i < np1; ++i) {
340 // Loop over 2 until a match is found.
341 bool match = false;
342 for (unsigned int j = 0; j < np2; ++j) {
343 if (fabs(xp2[j] - xp1[i]) < epsx && fabs(yp2[j] - yp1[i]) < epsy) {
344 match = true;
345 break;
346 }
347 const unsigned int jj = NextPoint(j, np2);
348 if (OnLine(xp2[j], yp2[j], xp2[jj], yp2[jj], xp1[i], yp1[i])) {
349 match = true;
350 break;
351 }
352 }
353 if (!match) return false;
354 }
355
356 // Compare all points of 2 with all points of 1.
357 for (unsigned int i = 0; i < np2; ++i) {
358 // Loop over 1 until a match is found.
359 bool match = false;
360 for (unsigned int j = 0; j < np1; ++j) {
361 if (fabs(xp2[i] - xp1[j]) < epsx && fabs(yp2[i] - yp1[j]) < epsy) {
362 match = true;
363 break;
364 }
365 const unsigned int jj = NextPoint(j, np1);
366 if (OnLine(xp1[j], yp1[j], xp1[jj], yp1[jj], xp2[i], yp2[i])) {
367 match = true;
368 break;
369 }
370 }
371 if (!match) return false;
372 }
373
374 // If we get this far, the curves are the same.
375 return true;
376}
377
378}
379
380namespace Garfield {
381
383
385
386Medium* ComponentNeBem3d::GetMedium(const double x, const double y,
387 const double z) {
388 if (!m_geometry) return nullptr;
389 return m_geometry->GetMedium(x, y, z, true);
390}
391
392void ComponentNeBem3d::ElectricField(const double x, const double y,
393 const double z, double& ex, double& ey,
394 double& ez, double& v, Medium*& m,
395 int& status) {
396 ex = ey = ez = v = 0.;
397 status = 0;
398 // Check if the requested point is inside a medium
399 m = GetMedium(x, y, z);
400 if (!m) {
401 status = -6;
402 } else if (!m->IsDriftable()) {
403 status = -5;
404 }
405
406 if (!m_ready) {
407 if (!Initialise()) {
408 std::cerr << m_className << "::ElectricField: Initialisation failed.\n";
409 status = -11;
410 return;
411 }
412 m_ready = true;
413 }
414
415 // Construct a point.
416 neBEM::Point3D point;
417 point.X = 0.01 * x;
418 point.Y = 0.01 * y;
419 point.Z = 0.01 * z;
420
421 // Compute the field.
422 neBEM::Vector3D field;
423 if (neBEM::neBEMField(&point, &v, &field) != 0) {
424 status = -10;
425 return;
426 }
427 ex = 0.01 * field.X;
428 ey = 0.01 * field.Y;
429 ez = 0.01 * field.Z;
430}
431
432void ComponentNeBem3d::ElectricField(const double x, const double y,
433 const double z, double& ex, double& ey,
434 double& ez, Medium*& m, int& status) {
435 double v = 0.;
436 ElectricField(x, y, z, ex, ey, ez, v, m, status);
437}
438
439bool ComponentNeBem3d::GetVoltageRange(double& vmin, double& vmax) {
440 // Voltage and other bc have to come from the solids
441 vmin = vmax = 0;
442 return true;
443}
444
445void ComponentNeBem3d::WeightingField(const double x, const double y,
446 const double z,
447 double& wx, double& wy, double& wz,
448 const std::string& label) {
449 wx = wy = wz = 0.;
450 if (m_wfields.count(label) == 0) return;
451 const int id = m_wfields[label];
452 neBEM::Point3D point;
453 point.X = 0.01 * x;
454 point.Y = 0.01 * y;
455 point.Z = 0.01 * z;
456 neBEM::Vector3D field;
457 if (neBEM::neBEMWeightingField(&point, &field, id) == DBL_MAX) {
458 std::cerr << m_className << "::WeightingField: Evaluation failed.\n";
459 return;
460 }
461 wx = 0.01 * field.X;
462 wy = 0.01 * field.Y;
463 wz = 0.01 * field.Z;
464}
465
466double ComponentNeBem3d::WeightingPotential(const double x, const double y,
467 const double z,
468 const std::string& label) {
469 if (m_wfields.count(label) == 0) return 0.;
470 const int id = m_wfields[label];
471 neBEM::Point3D point;
472 point.X = 0.01 * x;
473 point.Y = 0.01 * y;
474 point.Z = 0.01 * z;
475 neBEM::Vector3D field;
476 const double v = neBEM::neBEMWeightingField(&point, &field, id);
477 return v == DBL_MAX ? 0. : v;
478}
479
480void ComponentNeBem3d::AddPlaneX(const double x, const double v) {
481 if (m_ynplan[0] && m_ynplan[1]) {
482 std::cerr << m_className << "::AddPlaneX:\n"
483 << " Cannot have more than two planes at constant x.\n";
484 return;
485 }
486
487 if (m_ynplan[0]) {
488 m_ynplan[1] = true;
489 if (x < m_coplan[0]) {
490 m_coplan[1] = m_coplan[0];
491 m_vtplan[1] = m_vtplan[0];
492 m_coplan[0] = x;
493 m_vtplan[0] = v;
494 } else {
495 m_coplan[1] = x;
496 m_vtplan[1] = v;
497 }
498 } else {
499 m_ynplan[0] = true;
500 m_coplan[0] = x;
501 m_vtplan[0] = v;
502 }
503 m_ready = false;
504}
505
506void ComponentNeBem3d::AddPlaneY(const double y, const double v) {
507 if (m_ynplan[2] && m_ynplan[3]) {
508 std::cerr << m_className << "::AddPlaneY:\n"
509 << " Cannot have more than two planes at constant y.\n";
510 return;
511 }
512
513 if (m_ynplan[2]) {
514 m_ynplan[3] = true;
515 if (y < m_coplan[2]) {
516 m_coplan[3] = m_coplan[2];
517 m_vtplan[3] = m_vtplan[2];
518 m_coplan[2] = y;
519 m_vtplan[2] = v;
520 } else {
521 m_coplan[3] = y;
522 m_vtplan[3] = v;
523 }
524 } else {
525 m_ynplan[2] = true;
526 m_coplan[2] = y;
527 m_vtplan[2] = v;
528 }
529 m_ready = false;
530}
531
532void ComponentNeBem3d::AddPlaneZ(const double z, const double v) {
533 if (m_ynplan[4] && m_ynplan[5]) {
534 std::cerr << m_className << "::AddPlaneZ:\n"
535 << " Cannot have more than two planes at constant z.\n";
536 return;
537 }
538
539 if (m_ynplan[4]) {
540 m_ynplan[5] = true;
541 if (z < m_coplan[4]) {
542 m_coplan[5] = m_coplan[4];
543 m_vtplan[5] = m_vtplan[4];
544 m_coplan[4] = z;
545 m_vtplan[4] = v;
546 } else {
547 m_coplan[5] = z;
548 m_vtplan[5] = v;
549 }
550 } else {
551 m_ynplan[4] = true;
552 m_coplan[4] = z;
553 m_vtplan[4] = v;
554 }
555 m_ready = false;
556}
557
559 if (m_ynplan[0] && m_ynplan[1]) {
560 return 2;
561 } else if (m_ynplan[0] || m_ynplan[1]) {
562 return 1;
563 }
564 return 0;
565}
566
568 if (m_ynplan[2] && m_ynplan[3]) {
569 return 2;
570 } else if (m_ynplan[2] || m_ynplan[3]) {
571 return 1;
572 }
573 return 0;
574}
575
577 if (m_ynplan[4] && m_ynplan[5]) {
578 return 2;
579 } else if (m_ynplan[4] || m_ynplan[5]) {
580 return 1;
581 }
582 return 0;
583}
584
585bool ComponentNeBem3d::GetPlaneX(const unsigned int i, double& x,
586 double& v) const {
587 if (i >= 2 || (i == 1 && !m_ynplan[1])) {
588 std::cerr << m_className << "::GetPlaneX: Index out of range.\n";
589 return false;
590 }
591 x = m_coplan[i];
592 v = m_vtplan[i];
593 return true;
594}
595
596bool ComponentNeBem3d::GetPlaneY(const unsigned int i, double& y,
597 double& v) const {
598 if (i >= 2 || (i == 1 && !m_ynplan[3])) {
599 std::cerr << m_className << "::GetPlaneY: Index out of range.\n";
600 return false;
601 }
602 y = m_coplan[i + 2];
603 v = m_vtplan[i + 2];
604 return true;
605}
606
607bool ComponentNeBem3d::GetPlaneZ(const unsigned int i, double& z,
608 double& v) const {
609 if (i >= 2 || (i == 1 && !m_ynplan[5])) {
610 std::cerr << m_className << "::GetPlaneZ: Index out of range.\n";
611 return false;
612 }
613 z = m_coplan[i + 4];
614 v = m_vtplan[i + 4];
615 return true;
616}
617
618void ComponentNeBem3d::SetTargetElementSize(const double length) {
619
620 if (length < MinDist) {
621 std::cerr << m_className << "::SetTargetElementSize: Value must be > "
622 << MinDist << ".\n";
623 return;
624 }
625 m_targetElementSize = length;
626}
627
629 const unsigned int nmax) {
630
631 if (nmin == 0 || nmax == 0) {
632 std::cerr << m_className << "::SetMinMaxNumberOfElements:\n"
633 << " Values must be non-zero.\n";
634 return;
635 }
636 m_minNbElementsOnLength = std::min(nmin, nmax);
637 m_maxNbElementsOnLength = std::max(nmin, nmax);
638}
639
640void ComponentNeBem3d::SetPeriodicCopies(const unsigned int nx,
641 const unsigned int ny,
642 const unsigned int nz) {
643 m_nCopiesX = nx;
644 m_nCopiesY = ny;
645 m_nCopiesZ = nz;
646}
647
649 if (s < Small) {
650 std::cerr << m_className << "::SetPeriodicityX:\n"
651 << " Periodic length must be greater than zero.\n";
652 return;
653 }
654 m_periodicLength[0] = s;
655 m_periodic[0] = true;
656 m_mirrorPeriodic[0] = false;
658}
659
661 if (s < Small) {
662 std::cerr << m_className << "::SetPeriodicityY:\n"
663 << " Periodic length must be greater than zero.\n";
664 return;
665 }
666 m_periodicLength[1] = s;
667 m_periodic[1] = true;
668 m_mirrorPeriodic[1] = false;
670}
671
673 if (s < Small) {
674 std::cerr << m_className << "::SetPeriodicityZ:\n"
675 << " Periodic length must be greater than zero.\n";
676 return;
677 }
678 m_periodicLength[2] = s;
679 m_periodic[2] = true;
680 m_mirrorPeriodic[2] = false;
682}
683
685 if (s < Small) {
686 std::cerr << m_className << "::SetMirrorPeriodicityX:\n"
687 << " Periodic length must be greater than zero.\n";
688 return;
689 }
690 m_periodicLength[0] = s;
691 m_periodic[0] = false;
692 m_mirrorPeriodic[0] = true;
694}
695
697 if (s < Small) {
698 std::cerr << m_className << "::SetMirrorPeriodicityY:\n"
699 << " Periodic length must be greater than zero.\n";
700 return;
701 }
702 m_periodicLength[1] = s;
703 m_periodic[1] = false;
704 m_mirrorPeriodic[1] = true;
706}
707
709 if (s < Small) {
710 std::cerr << m_className << "::SetMirrorPeriodicityZ:\n"
711 << " Periodic length must be greater than zero.\n";
712 return;
713 }
714 m_periodicLength[2] = s;
715 m_periodic[2] = false;
716 m_mirrorPeriodic[2] = true;
718}
719
721 if (!m_periodic[0] && !m_mirrorPeriodic[0]) {
722 s = 0.;
723 return false;
724 }
725 s = m_periodicLength[0];
726 return true;
727}
728
730 if (!m_periodic[1] && !m_mirrorPeriodic[1]) {
731 s = 0.;
732 return false;
733 }
734 s = m_periodicLength[1];
735 return true;
736}
737
739 if (!m_periodic[2] && !m_mirrorPeriodic[2]) {
740 s = 0.;
741 return false;
742 }
743 s = m_periodicLength[2];
744 return true;
745}
746
748 // Reset the lists.
749 m_primitives.clear();
750 m_elements.clear();
751
752 if (!m_geometry) {
753 std::cerr << m_className << "::Initialise: Geometry not set.\n";
754 return false;
755 }
756 // Be sure we won't have intersections with the bounding box.
757 // TODO! Loop over the solids and call PLACYE, PLACHE, PLABXE, PLASPE, ...
758 // Loop over the solids.
759 std::map<int, Solid*> solids;
760 std::map<int, Solid::BoundaryCondition> bc;
761 std::map<int, double> volt;
762 std::map<int, double> eps;
763 std::map<int, double> charge;
764 const unsigned int nSolids = m_geometry->GetNumberOfSolids();
765 std::vector<Panel> panelsIn;
766 for (unsigned int i = 0; i < nSolids; ++i) {
767 Medium* medium = nullptr;
768 const auto solid = m_geometry->GetSolid(i, medium);
769 if (!solid) continue;
770 // Get the panels.
771 solid->SolidPanels(panelsIn);
772 // Get the boundary condition.
773 const auto id = solid->GetId();
774 solids[id] = solid;
775 bc[id] = solid->GetBoundaryConditionType();
776 volt[id] = solid->GetBoundaryPotential();
777 if (bc[id] == Solid::Unknown) {
778 std::cout << m_className << "::Initialise:\n"
779 << " Boundary conditions for solid " << id << " not set.\n";
780 if (medium && medium->IsConductor()) {
781 std::cout << " Assuming the panels to be grounded.\n";
782 bc[id] = Solid::Voltage;
783 volt[id] = 0.;
784 } else {
785 std::cout << " Assuming dielectric-dielectric interfaces.\n";
786 bc[id] = Solid::Dielectric;
787 }
788 }
789 charge[id] = solid->GetBoundaryChargeDensity();
790 if (!medium) {
791 eps[id] = 1.;
792 } else {
793 eps[id] = medium->GetDielectricConstant();
794 }
795 }
796 // Apply cuts.
797 // CALL CELSCT('APPLY')
798 // Reduce to basic periodic copy.
799 ShiftPanels(panelsIn);
800
801 // Find contact panels and split into primitives.
802
803 // *---------------------------------------------------------------------
804 // * PLABEM - Prepares panels for BEM applications: removes the contacts
805 // * and cuts polygons to rectangles and right-angle triangles.
806 // *---------------------------------------------------------------------
807
808 // Establish tolerances.
809 const double epsang = 1.e-6; // BEMEPA
810 const double epsxyz = 1.e-6; // BEMEPD
811 // CALL EPSSET('SET',EPSXYZ,EPSXYZ,EPSXYZ)
812
813 const unsigned int nIn = panelsIn.size();
814 if (m_debug) {
815 std::cout << m_className << "::Initialise: Retrieved "
816 << nIn << " panels from the solids.\n";
817 }
818 // Keep track of which panels have been processed.
819 std::vector<bool> mark(nIn, false);
820 // Count the number of interface panels that have been discarded.
821 unsigned int nTrivial = 0;
822 unsigned int nConflicting = 0;
823 unsigned int nNotImplemented = 0;
824 // Pick up panels which coincide potentially.
825 for (unsigned int i = 0; i < nIn; ++i) {
826 // Skip panels already done.
827 if (mark[i]) continue;
828 // Fetch panel parameters.
829 const double a1 = panelsIn[i].a;
830 const double b1 = panelsIn[i].b;
831 const double c1 = panelsIn[i].c;
832 const auto& xp1 = panelsIn[i].xv;
833 const auto& yp1 = panelsIn[i].yv;
834 const auto& zp1 = panelsIn[i].zv;
835 const unsigned int np1 = xp1.size();
836 // Establish its norm and offset.
837 const double d1 = a1 * xp1[0] + b1 * yp1[0] + c1 * zp1[0];
838 if (m_debug) {
839 std::cout << " Panel " << i << "\n Norm vector: "
840 << a1 << ", " << b1 << ", " << c1 << ", " << d1 << ".\n";
841 }
842 // Rotation matrix.
843 std::array<std::array<double, 3>, 3> rot;
844 if (fabs(c1) <= fabs(a1) && fabs(c1) <= fabs(b1)) {
845 // Rotation: removing C
846 rot[0][0] = b1 / sqrt(a1 * a1 + b1 * b1);
847 rot[0][1] = -a1 / sqrt(a1 * a1 + b1 * b1);
848 rot[0][2] = 0.;
849 } else if (fabs(b1) <= fabs(a1) && fabs(b1) <= fabs(c1)) {
850 // Rotation: removing B
851 rot[0][0] = c1 / sqrt(a1 * a1 + c1 * c1);
852 rot[0][1] = 0.;
853 rot[0][2] = -a1 / sqrt(a1 * a1 + c1 * c1);
854 } else {
855 // Rotation: removing A
856 rot[0][0] = 0.;
857 rot[0][1] = c1 / sqrt(b1 * b1 + c1 * c1);
858 rot[0][2] = -b1 / sqrt(b1 * b1 + c1 * c1);
859 }
860 rot[2][0] = a1;
861 rot[2][1] = b1;
862 rot[2][2] = c1;
863 rot[1][0] = rot[2][1] * rot[0][2] - rot[2][2] * rot[0][1];
864 rot[1][1] = rot[2][2] * rot[0][0] - rot[2][0] * rot[0][2];
865 rot[1][2] = rot[2][0] * rot[0][1] - rot[2][1] * rot[0][0];
866 // Rotate to the x, y plane.
867 std::vector<double> xp(np1, 0.);
868 std::vector<double> yp(np1, 0.);
869 std::vector<double> zp(np1, 0.);
870 double zm = 0.;
871 for (unsigned int k = 0; k < np1; ++k) {
872 xp[k] = rot[0][0] * xp1[k] + rot[0][1] * yp1[k] + rot[0][2] * zp1[k];
873 yp[k] = rot[1][0] * xp1[k] + rot[1][1] * yp1[k] + rot[1][2] * zp1[k];
874 zp[k] = rot[2][0] * xp1[k] + rot[2][1] * yp1[k] + rot[2][2] * zp1[k];
875 zm += zp[k];
876 }
877 zm /= np1;
878 // Store it.
879 std::vector<Panel> newPanels;
880 std::vector<int> vol1;
881 std::vector<int> vol2;
882 Panel panel1 = panelsIn[i];
883 panel1.xv = xp;
884 panel1.yv = yp;
885 panel1.zv = zp;
886 vol1.push_back(panel1.volume);
887 vol2.push_back(-1);
888 newPanels.push_back(std::move(panel1));
889 // Pick up all matching planes.
890 for (unsigned int j = i + 1; j < nIn; ++j) {
891 if (mark[j]) continue;
892 const double a2 = panelsIn[j].a;
893 const double b2 = panelsIn[j].b;
894 const double c2 = panelsIn[j].c;
895 const auto& xp2 = panelsIn[j].xv;
896 const auto& yp2 = panelsIn[j].yv;
897 const auto& zp2 = panelsIn[j].zv;
898 const unsigned int np2 = xp2.size();
899 // See whether this matches the first.
900 const double d2 = a2 * xp2[0] + b2 * yp2[0] + c2 * zp2[0];
901 // Inner product.
902 const double dot = a1 * a2 + b1 * b2 + c1 * c2;
903 // Offset between the two planes.
904 const double offset = d1 - d2 * dot;
905 if (fabs(fabs(dot) - 1.) > epsang || fabs(offset) > epsxyz) continue;
906 // Found a match.
907 mark[j] = true;
908 if (m_debug) std::cout << " Match with panel " << j << ".\n";
909 // Rotate this plane too.
910 xp.assign(np2, 0.);
911 yp.assign(np2, 0.);
912 zp.assign(np2, 0.);
913 zm = 0.;
914 for (unsigned int k = 0; k < np2; ++k) {
915 xp[k] = rot[0][0] * xp2[k] + rot[0][1] * yp2[k] + rot[0][2] * zp2[k];
916 yp[k] = rot[1][0] * xp2[k] + rot[1][1] * yp2[k] + rot[1][2] * zp2[k];
917 zp[k] = rot[2][0] * xp2[k] + rot[2][1] * yp2[k] + rot[2][2] * zp2[k];
918 zm += zp[k];
919 }
920 zm /= np2;
921 // Store it.
922 Panel panel2 = panelsIn[j];
923 panel2.xv = xp;
924 panel2.yv = yp;
925 panel2.zv = zp;
926 vol1.push_back(panel2.volume);
927 vol2.push_back(-1);
928 newPanels.push_back(std::move(panel2));
929 }
930 std::vector<bool> obsolete(newPanels.size(), false);
931 // Cut them as long as needed till no contacts remain.
932 unsigned int jmin = 0;
933 bool change = true;
934 while (change) {
935 change = false;
936 const unsigned int n = newPanels.size();
937 for (unsigned int j = 0; j < n; ++j) {
938 if (obsolete[j] || j < jmin) continue;
939 if (vol1[j] >= 0 && vol2[j] >= 0) continue;
940 const auto& panelj = newPanels[j];
941 for (unsigned int k = j + 1; k < n; ++k) {
942 if (obsolete[k]) continue;
943 if (vol1[k] >= 0 && vol2[k] >= 0) continue;
944 const auto& panelk = newPanels[k];
945 if (m_debug) std::cout << " Cutting " << j << ", " << k << ".\n";
946 // Separate contact and non-contact areas.
947 std::vector<Panel> panelsOut;
948 std::vector<int> itypo;
949 EliminateOverlaps(panelj, panelk, panelsOut, itypo);
950 const unsigned int nOut = panelsOut.size();
951 if (nOut == 2) {
952 // TODO: retrieve epsx, epsy from overlap finding?
953 const double epsx = epsxyz;
954 const double epsy = epsxyz;
955 // If there are just 2 panels, see whether there is a new one.
956 const bool equal1 = Equal(panelj, panelsOut[0], epsx, epsy);
957 const bool equal2 = Equal(panelj, panelsOut[1], epsx, epsy);
958 const bool equal3 = Equal(panelk, panelsOut[0], epsx, epsy);
959 const bool equal4 = Equal(panelk, panelsOut[1], epsx, epsy);
960 if ((equal1 || equal3) && (equal2 || equal4)) {
961 if (m_debug) {
962 std::cout << " Original and new panels are identical.\n";
963 }
964 } else {
965 change = true;
966 }
967 } else {
968 change = true;
969 }
970 if (m_debug) std::cout << " Change flag: " << change << ".\n";
971 // If there is no change, keep the two panels and proceed.
972 if (!change) continue;
973 // Flag the existing panels as inactive.
974 obsolete[j] = true;
975 obsolete[k] = true;
976
977 // Add the new panels.
978 for (unsigned int l = 0; l < nOut; ++l) {
979 if (itypo[l] == 1) {
980 vol1.push_back(std::max(vol1[j], vol2[j]));
981 vol2.push_back(-1);
982 } else if (itypo[l] == 2) {
983 vol1.push_back(std::max(vol1[k], vol2[k]));
984 vol2.push_back(-1);
985 } else {
986 vol1.push_back(std::max(vol1[j], vol2[j]));
987 vol2.push_back(std::max(vol1[k], vol2[k]));
988 }
989 newPanels.push_back(std::move(panelsOut[l]));
990 obsolete.push_back(false);
991 }
992 jmin = j + 1;
993 // Restart the loops.
994 break;
995 }
996 if (change) break;
997 }
998 }
999 // And rotate the panels back in place.
1000 const unsigned int nNew = newPanels.size();
1001 for (unsigned int j = 0; j < nNew; ++j) {
1002 if (obsolete[j]) continue;
1003 // Examine the boundary conditions.
1004 int interfaceType = 0;
1005 double potential = 0.;
1006 double lambda = 0.;
1007 double chargeDensity = 0.;
1008 if (m_debug) {
1009 std::cout << " Volume 1: " << vol1[j]
1010 << ". Volume 2: " << vol2[j] << ".\n";
1011 }
1012 if (vol1[j] < 0 && vol2[j] < 0) {
1013 // Shouldn't happen.
1014 continue;
1015 } else if (vol1[j] < 0 || vol2[j] < 0) {
1016 // Interface between a solid and vacuum/background.
1017 const auto vol = vol1[j] < 0 ? vol2[j] : vol1[j];
1018 interfaceType = InterfaceType(bc[vol]);
1019 if (bc[vol] == Solid::Dielectric ||
1020 bc[vol] == Solid::DielectricCharge) {
1021 if (fabs(eps[vol] - 1.) < 1.e-6) {
1022 // Same epsilon on both sides. Skip.
1023 interfaceType = 0;
1024 } else {
1025 lambda = (eps[vol] - 1.) / (eps[vol] + 1.);
1026 }
1027 } else if (bc[vol] == Solid::Voltage) {
1028 potential = volt[vol];
1029 }
1030 if (bc[vol] == Solid::Charge || bc[vol] == Solid::DielectricCharge) {
1031 chargeDensity = charge[vol];
1032 }
1033 } else {
1034 const auto bc1 = bc[vol1[j]];
1035 const auto bc2 = bc[vol2[j]];
1036 if (bc1 == Solid::Voltage || bc1 == Solid::Charge ||
1037 bc1 == Solid::Float) {
1038 interfaceType = InterfaceType(bc1);
1039 // First volume is a conductor. Other volume must be a dielectric.
1040 if (bc2 == Solid::Dielectric || bc2 == Solid::DielectricCharge) {
1041 if (bc1 == Solid::Voltage) {
1042 potential = volt[vol1[j]];
1043 } else if (bc1 == Solid::Charge) {
1044 chargeDensity = charge[vol1[j]];
1045 }
1046 } else {
1047 interfaceType = -1;
1048 }
1049 if (bc1 == Solid::Voltage && bc2 == Solid::Voltage) {
1050 const double v1 = volt[vol1[j]];
1051 const double v2 = volt[vol2[j]];
1052 if (fabs(v1 - v2) < 1.e-6 * (1. + fabs(v1) + fabs(v2))) {
1053 interfaceType = 0;
1054 }
1055 }
1056 } else if (bc1 == Solid::Dielectric ||
1057 bc1 == Solid::DielectricCharge) {
1058 interfaceType = InterfaceType(bc2);
1059 // First volume is a dielectric.
1060 if (bc2 == Solid::Voltage) {
1061 potential = volt[vol2[j]];
1062 } else if (bc2 == Solid::Charge) {
1063 chargeDensity = charge[vol2[j]];
1064 } else if (bc2 == Solid::Dielectric ||
1065 bc2 == Solid::DielectricCharge) {
1066 const double eps1 = eps[vol1[j]];
1067 const double eps2 = eps[vol2[j]];
1068 if (fabs(eps1 - eps2) < 1.e-6 * (1. + fabs(eps1) + fabs(eps2))) {
1069 // Same epsilon. Skip.
1070 interfaceType = 0;
1071 } else {
1072 lambda = (eps1 - eps2) / (eps1 + eps2);
1073 }
1074 }
1075 }
1076 }
1077 if (m_debug) {
1078 if (interfaceType < 0) {
1079 std::cout << " Conflicting boundary conditions. Skip.\n";
1080 } else if (interfaceType < 1) {
1081 std::cout << " Trivial interface. Skip.\n";
1082 } else if (interfaceType > 5) {
1083 std::cout << " Interface type " << interfaceType
1084 << " is not implemented. Skip.\n";
1085 }
1086 }
1087 if (interfaceType < 0) {
1088 ++nConflicting;
1089 } else if (interfaceType < 1) {
1090 ++nTrivial;
1091 } else if (interfaceType > 5) {
1092 ++nNotImplemented;
1093 }
1094 if (interfaceType < 1 || interfaceType > 5) continue;
1095
1096 std::vector<Panel> panelsOut;
1097 // Reduce to rectangles and right-angle triangles.
1098 if (m_debug) std::cout << " Creating primitives.\n";
1099 MakePrimitives(newPanels[j], panelsOut);
1100 // Loop over the rectangles and triangles.
1101 for (auto& panel : panelsOut) {
1102 const auto& up = panel.xv;
1103 const auto& vp = panel.yv;
1104 const auto& wp = panel.zv;
1105 const unsigned int np = up.size();
1106 // Rotate.
1107 xp.assign(np, 0.);
1108 yp.assign(np, 0.);
1109 zp.assign(np, 0.);
1110 for (unsigned int k = 0; k < np; ++k) {
1111 xp[k] = rot[0][0] * up[k] + rot[1][0] * vp[k] + rot[2][0] * wp[k];
1112 yp[k] = rot[0][1] * up[k] + rot[1][1] * vp[k] + rot[2][1] * wp[k];
1113 zp[k] = rot[0][2] * up[k] + rot[1][2] * vp[k] + rot[2][2] * wp[k];
1114 }
1115 Primitive primitive;
1116 primitive.a = panel.a;
1117 primitive.b = panel.b;
1118 primitive.c = panel.c;
1119 primitive.xv = xp;
1120 primitive.yv = yp;
1121 primitive.zv = zp;
1122 primitive.v = potential;
1123 primitive.q = chargeDensity;
1124 primitive.lambda = lambda;
1125 primitive.interface = interfaceType;
1126 // Set the requested discretization level (target element size).
1127 primitive.elementSize = -1.;
1128 if (solids.find(vol1[j]) != solids.end()) {
1129 const auto solid = solids[vol1[j]];
1130 if (solid) {
1131 primitive.elementSize = solid->GetDiscretisationLevel(panel);
1132 }
1133 }
1134 primitive.vol1 = vol1[j];
1135 primitive.vol2 = vol2[j];
1136 m_primitives.push_back(std::move(primitive));
1137 }
1138 }
1139 }
1140
1141 // Add the wires.
1142 for (unsigned int i = 0; i < nSolids; ++i) {
1143 const auto solid = m_geometry->GetSolid(i);
1144 if (!solid) continue;
1145 if (!solid->IsWire()) continue;
1146 double x0 = 0., y0 = 0., z0 = 0.;
1147 solid->GetCentre(x0, y0, z0);
1148 double dx = 0., dy = 0., dz = 1.;
1149 solid->GetDirection(dx, dy, dz);
1150 const double dnorm = sqrt(dx * dx + dy * dy + dz * dz);
1151 if (dnorm < Small) {
1152 std::cerr << m_className << "::Initialise:\n"
1153 << " Wire has zero norm direction vector; skipped.\n";
1154 continue;
1155 }
1156 dx /= dnorm;
1157 dy /= dnorm;
1158 dz /= dnorm;
1159 const double h = solid->GetHalfLengthZ();
1160 Primitive primitive;
1161 primitive.a = solid->GetRadius();
1162 primitive.b = 0.;
1163 primitive.c = 0.;
1164 primitive.xv = {x0 - h * dx, x0 + h * dx};
1165 primitive.yv = {y0 - h * dy, y0 + h * dy};
1166 primitive.zv = {z0 - h * dz, z0 + h * dz};
1167 primitive.v = solid->GetBoundaryPotential();
1168 primitive.q = solid->GetBoundaryChargeDensity();
1169 primitive.lambda = 0.;
1170 primitive.interface = InterfaceType(solid->GetBoundaryConditionType());
1171 // Set the requested discretization level (target element size).
1172 Panel panel;
1173 primitive.elementSize = solid->GetDiscretisationLevel(panel);
1174 primitive.vol1 = solid->GetId();
1175 primitive.vol2 = -1;
1176 m_primitives.push_back(std::move(primitive));
1177 }
1178 // Print a warning if we have discarded some panels during the process.
1179 if (nTrivial > 0 || nConflicting > 0 || nNotImplemented > 0) {
1180 std::cerr << m_className << "::Initialise:\n";
1181 if (nConflicting > 0) {
1182 std::cerr << " Skipped " << nConflicting
1183 << " panels with conflicting boundary conditions.\n";
1184 }
1185 if (nNotImplemented > 0) {
1186 std::cerr << " Skipped " << nNotImplemented
1187 << " panels with not yet available boundary conditions.\n";
1188 }
1189 if (nTrivial > 0) {
1190 std::cerr << " Skipped " << nTrivial
1191 << " panels with trivial boundary conditions.\n";
1192 }
1193 }
1194 if (m_debug) {
1195 std::cout << m_className << "::Initialise:\n"
1196 << " Created " << m_primitives.size() << " primitives.\n";
1197 }
1198
1199 // Discretize the primitives.
1200 for (const auto& primitive : m_primitives) {
1201 const auto nVertices = primitive.xv.size();
1202 if (nVertices < 2 || nVertices > 4) continue;
1203 std::vector<Element> elements;
1204 // Get the target element size.
1205 double targetSize = primitive.elementSize;
1206 if (targetSize < MinDist) targetSize = m_targetElementSize;
1207 if (nVertices == 2) {
1208 DiscretizeWire(primitive, targetSize, elements);
1209 } else if (nVertices == 3) {
1210 DiscretizeTriangle(primitive, targetSize, elements);
1211 } else if (nVertices == 4) {
1212 DiscretizeRectangle(primitive, targetSize, elements);
1213 }
1214 for (auto& element: elements) {
1215 element.interface = primitive.interface;
1216 element.lambda = primitive.lambda;
1217 element.bc = primitive.v;
1218 element.assigned = primitive.q;
1219 element.solution = 0.;
1220 }
1221 m_elements.insert(m_elements.end(),
1222 std::make_move_iterator(elements.begin()),
1223 std::make_move_iterator(elements.end()));
1224 }
1225 neBEM::NbThreads = m_nThreads;
1226 if (neBEM::neBEMInitialize() != 0) {
1227 std::cerr << m_className << "::Initialise:\n"
1228 << " Initialising neBEM failed.\n";
1229 return false;
1230 }
1231 gComponentNeBem3d = this;
1232 // Set the user options.
1233 neBEM::MinNbElementsOnLength = m_minNbElementsOnLength;
1234 neBEM::MaxNbElementsOnLength = m_maxNbElementsOnLength;
1235 neBEM::ElementLengthRqstd = m_targetElementSize * 0.01;
1236
1237 // New model / reuse existing model flag.
1238 // TODO!
1239 neBEM::NewModel = 1;
1240 neBEM::NewMesh = 1;
1241 neBEM::NewBC = 1;
1242 neBEM::NewPP = 1;
1243 // Pass debug level.
1244 neBEM::DebugLevel = m_debug ? 101 : 0;
1245 // Store inverted matrix or not.
1246 // TODO!
1247 neBEM::OptStoreInvMatrix = 0;
1248 neBEM::OptFormattedFile = 0;
1249 // Matrix inversion method (LU or SVD).
1250 if (m_inversion == Inversion::LU) {
1251 neBEM::OptInvMatProc = 0;
1252 } else {
1253 neBEM::OptInvMatProc = 1;
1254 }
1255 // Delete existing weighting fields (if any).
1256 if (neBEM::WtFieldChDen != NULL) {
1257 neBEM::neBEMDeleteAllWeightingFields();
1258 }
1259 // Transfer the geometry.
1260 if (neBEM::neBEMReadGeometry() != 0) {
1261 std::cerr << m_className << "::Initialise:\n"
1262 << " Transferring the geometry to neBEM failed.\n";
1263 return false;
1264 }
1265
1266 // Discretization.
1267 int** elementNbs = neBEM::imatrix(1, neBEM::NbPrimitives, 1, 2);
1268 for (int i = 1; i <= neBEM::NbPrimitives; ++i) {
1269 const int vol1 = neBEM::VolRef1[i];
1270 double size1 = -1.;
1271 if (solids.find(vol1) != solids.end()) {
1272 const auto solid = solids[vol1];
1273 if (solid) {
1274 Panel panel;
1275 panel.a = neBEM::XNorm[i];
1276 panel.b = neBEM::YNorm[i];
1277 panel.c = neBEM::ZNorm[i];
1278 const int nv = neBEM::NbVertices[i];
1279 panel.xv.resize(nv);
1280 panel.yv.resize(nv);
1281 panel.zv.resize(nv);
1282 for (int j = 0; j < nv; ++j) {
1283 panel.xv[j] = neBEM::XVertex[i][j];
1284 panel.yv[j] = neBEM::YVertex[i][j];
1285 panel.zv[j] = neBEM::ZVertex[i][j];
1286 }
1287 size1 = solid->GetDiscretisationLevel(panel);
1288 }
1289 }
1290 int vol2 = neBEM::VolRef2[i];
1291 double size2 = -1.;
1292 if (solids.find(vol2) != solids.end()) {
1293 const auto solid = solids[vol2];
1294 if (solid) {
1295 Panel panel;
1296 panel.a = -neBEM::XNorm[i];
1297 panel.b = -neBEM::YNorm[i];
1298 panel.c = -neBEM::ZNorm[i];
1299 const int nv = neBEM::NbVertices[i];
1300 panel.xv.resize(nv);
1301 panel.yv.resize(nv);
1302 panel.zv.resize(nv);
1303 for (int j = 0; j < nv; ++j) {
1304 panel.xv[j] = neBEM::XVertex[i][j];
1305 panel.yv[j] = neBEM::YVertex[i][j];
1306 panel.zv[j] = neBEM::ZVertex[i][j];
1307 }
1308 size2 = solid->GetDiscretisationLevel(panel);
1309 }
1310 }
1311 double size = m_targetElementSize;
1312 if (size1 > 0. && size2 > 0.) {
1313 size = std::min(size1, size2);
1314 } else if (size1 > 0.) {
1315 size = size1;
1316 } else if (size2 > 0.) {
1317 size = size2;
1318 }
1319 // Convert from cm to m.
1320 size *= 0.01;
1321
1322 // Work out the element dimensions.
1323 const double dx1 = neBEM::XVertex[i][0] - neBEM::XVertex[i][1];
1324 const double dy1 = neBEM::YVertex[i][0] - neBEM::YVertex[i][1];
1325 const double dz1 = neBEM::ZVertex[i][0] - neBEM::ZVertex[i][1];
1326 const double l1 = dx1 * dx1 + dy1 * dy1 + dz1 * dz1;
1327 int nb1 = (int)(sqrt(l1) / size);
1328 // Truncate to desired range.
1329 if (nb1 < neBEM::MinNbElementsOnLength) {
1330 nb1 = neBEM::MinNbElementsOnLength;
1331 } else if (nb1 > neBEM::MaxNbElementsOnLength) {
1332 nb1 = neBEM::MaxNbElementsOnLength;
1333 }
1334 int nb2 = 0;
1335 if (neBEM::NbVertices[i] > 2) {
1336 const double dx2 = neBEM::XVertex[i][2] - neBEM::XVertex[i][1];
1337 const double dy2 = neBEM::YVertex[i][2] - neBEM::YVertex[i][1];
1338 const double dz2 = neBEM::ZVertex[i][2] - neBEM::ZVertex[i][1];
1339 const double l2 = dx2 * dx2 + dy2 * dy2 + dz2 * dz2;
1340 nb2 = (int)(sqrt(l2) / size);
1341 if (nb2 < neBEM::MinNbElementsOnLength) {
1342 nb2 = neBEM::MinNbElementsOnLength;
1343 } else if (nb2 > neBEM::MaxNbElementsOnLength) {
1344 nb2 = neBEM::MaxNbElementsOnLength;
1345 }
1346 }
1347 elementNbs[i][1] = nb1;
1348 elementNbs[i][2] = nb2;
1349 }
1350
1351 if (neBEM::neBEMDiscretize(elementNbs) != 0) {
1352 std::cerr << m_className << "::Initialise: Discretization failed.\n";
1353 neBEM::free_imatrix(elementNbs, 1, neBEM::NbPrimitives, 1, 2);
1354 return false;
1355 }
1356 neBEM::free_imatrix(elementNbs, 1, neBEM::NbPrimitives, 1, 2);
1357 if (neBEM::neBEMBoundaryConditions() != 0) {
1358 std::cerr << m_className << "::Initialise:\n"
1359 << " Setting the boundary conditions failed.\n";
1360 return false;
1361 }
1362 if (neBEM::neBEMSolve() != 0) {
1363 std::cerr << m_className << "::Initialise: Solution failed.\n";
1364 return false;
1365 }
1366 // Now the weighting fields.
1367 std::set<std::string> labels;
1368 for (unsigned int i = 0; i < nSolids; ++i) {
1369 const auto solid = m_geometry->GetSolid(i);
1370 if (!solid) continue;
1371 const std::string label = solid->GetLabel();
1372 if (!label.empty()) labels.insert(label);
1373 }
1374 for (const auto& label : labels) {
1375 std::vector<int> primitives;
1376 for (unsigned int i = 0; i < nSolids; ++i) {
1377 const auto solid = m_geometry->GetSolid(i);
1378 if (!solid) continue;
1379 if (solid->GetLabel() != label) continue;
1380 const int id = solid->GetId();
1381 // Add the primitives associated to this solid to the list.
1382 for (int j = 1; j <= neBEM::NbPrimitives; ++j) {
1383 if (neBEM::VolRef1[j] == id || neBEM::VolRef2[j] == id) {
1384 primitives.push_back(j);
1385 }
1386 }
1387 }
1388 // Request the weighting field for this list of primitives.
1389 const int np = primitives.size();
1390 const int id = neBEM::neBEMPrepareWeightingField(np, primitives.data());
1391 if (id < 0) {
1392 std::cerr << m_className << "::Initialise:\n"
1393 << " Weighting field calculation for readout group \""
1394 << label << "\" failed.\n";
1395 continue;
1396 } else {
1397 std::cout << m_className << "::Initialise:\n"
1398 << " Prepared weighting field for readout group \""
1399 << label << "\".\n";
1400 m_wfields[label] = id;
1401 }
1402 }
1403 // TODO! Not sure if we should call this here.
1404 // neBEM::neBEMEnd();
1405 m_ready = true;
1406 return true;
1407}
1408
1409void ComponentNeBem3d::ShiftPanels(std::vector<Panel>& panels) const {
1410
1411 // *---------------------------------------------------------------------
1412 // * BEMBAS - Reduces panels to the basic period.
1413 // *---------------------------------------------------------------------
1414 const bool perx = m_periodic[0] || m_mirrorPeriodic[0];
1415 const bool pery = m_periodic[1] || m_mirrorPeriodic[1];
1416 const bool perz = m_periodic[2] || m_mirrorPeriodic[2];
1417 // Nothing to do if there is no periodicity.
1418 if (!perx && !pery && !perz) return;
1419
1420 // Loop over the panels.
1421 for (auto& panel : panels) {
1422 const auto nv = panel.xv.size();
1423 if (nv == 0) continue;
1424 // Determine the centre of gravity.
1425 double xc = std::accumulate(panel.xv.begin(), panel.xv.end(), 0.);
1426 double yc = std::accumulate(panel.yv.begin(), panel.yv.end(), 0.);
1427 double zc = std::accumulate(panel.zv.begin(), panel.zv.end(), 0.);
1428 xc /= nv;
1429 yc /= nv;
1430 zc /= nv;
1431 // Any change ?
1432 constexpr double eps = 1.e-6;
1433 double rx = 0.;
1434 int nx = 0;
1435 if (perx && m_periodicLength[0] > 0.) {
1436 rx = xc / m_periodicLength[0];
1437 nx = std::round(rx);
1438 if (std::abs(rx - nx - 0.5) < eps) ++nx;
1439 }
1440 double ry = 0.;
1441 int ny = 0;
1442 if (pery && m_periodicLength[1] > 0.) {
1443 ry = yc / m_periodicLength[1];
1444 ny = std::round(ry);
1445 if (std::abs(ry - ny - 0.5) < eps) ++ny;
1446 }
1447 double rz = 0.;
1448 int nz = 0;
1449 if (perz && m_periodicLength[2] > 0.) {
1450 rz = zc / m_periodicLength[2];
1451 nz = std::round(rz);
1452 if (std::abs(rz - nz - 0.5) < eps) ++nz;
1453 }
1454 // Skip if there is no shift.
1455 if (nx == 0 && ny == 0 && nz == 0) continue;
1456 // Shift for x-periodicity.
1457 if (nx != 0) {
1458 const double shift = nx * m_periodicLength[0];
1459 for (auto& x : panel.xv) x -= shift;
1460 }
1461 // Shift for y-periodicity.
1462 if (ny != 0) {
1463 const double shift = ny * m_periodicLength[1];
1464 for (auto& y : panel.yv) y -= shift;
1465 }
1466 // Shift for z-periodicity.
1467 if (nz != 0) {
1468 const double shift = nz * m_periodicLength[2];
1469 for (auto& z : panel.zv) z -= shift;
1470 }
1471 }
1472}
1473
1474int ComponentNeBem3d::InterfaceType(const Solid::BoundaryCondition bc) const {
1475
1476 // 0: To be skipped
1477 // 1: Conductor-dielectric
1478 // 2: Conductor with known charge
1479 // 3: Conductor at floating potential
1480 // 4: Dielectric-dielectric
1481 // 5: Dielectric with given surface charge
1482 //
1483 // Check dielectric-dielectric formulation in
1484 // Jaydeep P. Bardhan,
1485 // Numerical solution of boundary-integral equations for molecular
1486 // electrostatics, J. Chem. Phys. 130, 094102 (2009)
1487 // https://doi.org/10.1063/1.3080769
1488
1489 switch (bc) {
1490 case Solid::Voltage:
1491 return 1;
1492 break;
1493 case Solid::Charge:
1494 return 2;
1495 break;
1496 case Solid::Float:
1497 return 3;
1498 break;
1499 case Solid::Dielectric:
1500 return 4;
1501 break;
1503 return 5;
1504 break;
1506 // Symmetry boundary, E parallel.
1507 return 6;
1508 break;
1510 // Symmetry boundary, E perpendicular.
1511 return 7;
1512 break;
1513 default:
1514 break;
1515 }
1516 return 0;
1517}
1518
1519bool ComponentNeBem3d::DiscretizeWire(const Primitive& primitive,
1520 const double targetSize, std::vector<Element>& elements) const {
1521
1522 const double dx = primitive.xv[1] - primitive.xv[0];
1523 const double dy = primitive.yv[1] - primitive.yv[0];
1524 const double dz = primitive.zv[1] - primitive.zv[0];
1525 const double lw = sqrt(dx * dx + dy * dy + dz * dz);
1526 // Determine the number of segments along the wire.
1527 unsigned int nSegments = NbOfSegments(lw, targetSize);
1528 const double elementSize = lw / nSegments;
1529
1530 // Determine the direction cosines.
1531 // The z axis of the local coordinate system is along the wire.
1532 const std::array<double, 3> zu = {dx / lw, dy / lw, dz / lw};
1533 // Make the x axis orthogonal in the two largest components.
1534 std::array<double, 3> xu;
1535 if (fabs(zu[0]) >= fabs(zu[2]) && fabs(zu[1]) >= fabs(zu[2])) {
1536 xu = {-zu[1], zu[0], 0.};
1537 } else if (fabs(zu[0]) >= fabs(zu[1]) && fabs(zu[2]) >= fabs(zu[1])) {
1538 xu = {-zu[2], 0., zu[0]};
1539 } else {
1540 xu = {0., zu[2], -zu[1]};
1541 }
1542 xu = UnitVector(xu);
1543 // The y axis is given by the vectorial product of z and x.
1544 const std::array<double, 3> yu = UnitVector(CrossProduct(zu, xu));
1545 const std::array<std::array<double, 3>, 3> dcos = {xu, yu, zu};
1546
1547 const double xincr = dx / nSegments;
1548 const double yincr = dy / nSegments;
1549 const double zincr = dz / nSegments;
1550
1551 // TODO!
1552 const double radius = 1.;
1553 const double dA = TwoPi * radius * elementSize;
1554 for (unsigned int i = 0; i < nSegments; ++i) {
1555 const double x0 = primitive.xv[0] + i * xincr;
1556 const double y0 = primitive.yv[0] + i * yincr;
1557 const double z0 = primitive.zv[0] + i * zincr;
1558 Element element;
1559 element.origin = {x0 + 0.5 * xincr, y0 + 0.5 * yincr, z0 + 0.5 * zincr};
1560 element.xv = {x0, x0 + xincr};
1561 element.yv = {y0, y0 + yincr};
1562 element.zv = {z0, z0 + zincr};
1563 element.lx = radius;
1564 element.lz = elementSize;
1565 element.dA = dA;
1566 element.dcos = dcos;
1567 // Modify to be on surface?
1568 element.collocationPoint = element.origin;
1569 elements.push_back(std::move(element));
1570 }
1571 return true;
1572}
1573
1574bool ComponentNeBem3d::DiscretizeTriangle(const Primitive& primitive,
1575 const double targetSize, std::vector<Element>& elements) const {
1576
1577 // Origin of the local coordinate system is at the right angle corner.
1578 std::array<double, 3> corner = {primitive.xv[1], primitive.yv[1],
1579 primitive.zv[1]};
1580 // Determine the direction cosines.
1581 const double dx1 = primitive.xv[0] - primitive.xv[1];
1582 const double dy1 = primitive.yv[0] - primitive.yv[1];
1583 const double dz1 = primitive.zv[0] - primitive.zv[1];
1584 const double dx2 = primitive.xv[2] - primitive.xv[1];
1585 const double dy2 = primitive.yv[2] - primitive.yv[1];
1586 const double dz2 = primitive.zv[2] - primitive.zv[1];
1587 // Normal vector of the primitive.
1588 std::array<double, 3> nu = {primitive.a, primitive.b, primitive.c};
1589 nu = UnitVector(nu);
1590 // int flagDC = 1;
1591 // We begin with trial 1: one of the possible orientations.
1592 double lx = sqrt(dx1 * dx1 + dy1 * dy1 + dz1 * dz1);
1593 double lz = sqrt(dx2 * dx2 + dy2 * dy2 + dz2 * dz2);
1594 std::array<double, 3> xu = {dx1 / lx, dy1 / lx, dz1 / lx};
1595 std::array<double, 3> zu = {dx2 / lz, dy2 / lz, dz2 / lz};
1596 std::array<double, 3> yu = CrossProduct(zu, xu);
1597 constexpr double tol = 1.e-3;
1598 if ((fabs(yu[0] - nu[0]) > tol) || (fabs(yu[1] - nu[1]) > tol) ||
1599 (fabs(yu[2] - nu[2]) > tol)) {
1600 // flagDC = 2;
1601 // Try the other orientation.
1602 std::swap(lx, lz);
1603 xu = {dx2 / lx, dy2 / lx, dz2 / lx};
1604 zu = {dx1 / lz, dy1 / lz, dz1 / lz};
1605 yu = CrossProduct(zu, xu);
1606 if ((fabs(yu[0] - nu[0]) > tol) || (fabs(yu[1] - nu[1]) > tol) ||
1607 (fabs(yu[2] - nu[2]) > tol)) {
1608 // No other possibility, search failed.
1609 std::cerr << m_className << "::DiscretizeTriangle:\n"
1610 << " Could not establish direction vectors.\n";
1611 return false;
1612 }
1613 }
1614 std::array<std::array<double, 3>, 3> dcos = {xu, yu, zu};
1615
1616 // TODO!
1617 /*
1618 double eps1 = primitive.eps1;
1619 double eps2 = primitive.eps2;
1620 if (flagDC == 1) std::swap(eps1, eps2);
1621 */
1622
1623 // Determine the number of elements.
1624 unsigned int nx = NbOfSegments(lx, targetSize);
1625 unsigned int nz = NbOfSegments(lz, targetSize);
1626 double elementSizeX = lx / nx;
1627 double elementSizeZ = lz / nz;
1628
1629 // Analyse the aspect ratio.
1630 double ar = elementSizeX / elementSizeZ;
1631 if (ar > 10.) {
1632 // Reduce nz.
1633 elementSizeZ = 0.1 * elementSizeX;
1634 nz = std::max(static_cast<int>(lz / elementSizeZ), 1);
1635 elementSizeZ = lz / nz;
1636 }
1637 if (ar < 0.1) {
1638 // Reduce nx.
1639 elementSizeX = 0.1 * elementSizeZ;
1640 nx = std::max(static_cast<int>(lx / elementSizeX), 1);
1641 elementSizeX = lx / nx;
1642 }
1643 const double dxdz = lx / lz;
1644 for (unsigned int k = 0; k < nz; ++k) {
1645 // Consider the k-th row.
1646 const double zlo = k * elementSizeZ;
1647 const double zhi = (k + 1) * elementSizeZ;
1648 double xlo = (lz - zlo) * dxdz;
1649 double xhi = (lz - zhi) * dxdz;
1650 // Triangular element on the k-th row.
1651 Element triangle;
1652 triangle.origin = LocalToGlobal(xhi, 0., zlo, dcos, corner);
1653 // Assign element values.
1654 const double triangleSizeX = xlo - xhi;
1655 triangle.lx = triangleSizeX;
1656 triangle.lz = elementSizeZ;
1657 triangle.dA = 0.5 * triangleSizeX * elementSizeZ;
1658 triangle.dcos = dcos;
1659 // Calculate the element vertices.
1660 const double xv0 = triangle.origin[0];
1661 const double yv0 = triangle.origin[1];
1662 const double zv0 = triangle.origin[2];
1663 const double xv1 = xv0 + triangleSizeX * xu[0];
1664 const double yv1 = yv0 + triangleSizeX * xu[1];
1665 const double zv1 = zv0 + triangleSizeX * xu[2];
1666 const double xv2 = xv0 + elementSizeZ * zu[0];
1667 const double yv2 = yv0 + elementSizeZ * zu[1];
1668 const double zv2 = zv0 + elementSizeZ * zu[2];
1669 // Assign vertices of the element.
1670 triangle.xv = {xv0, xv1, xv2};
1671 triangle.yv = {yv0, yv1, yv2};
1672 triangle.zv = {zv0, zv1, zv2};
1673
1674 // Boundary condition is applied at the barycentre, not at the origin
1675 // of the element coordinate system which is at the right-angled corner.
1676 const double xb = triangleSizeX / 3.;
1677 const double yb = 0.;
1678 const double zb = elementSizeZ / 3.;
1679 triangle.collocationPoint = LocalToGlobal(xb, yb, zb, dcos, triangle.origin);
1680 elements.push_back(std::move(triangle));
1681 // No rectangular element on the last row.
1682 if (k == nz - 1) continue;
1683
1684 // Determine number and size in x of the rectangular elements.
1685 const int nRect = xhi <= elementSizeX ? 1 : (int)(xhi / elementSizeX);
1686 const double rectSizeX = xhi / nRect;
1687 const double zc = 0.5 * (zlo + zhi);
1688 for (int i = 0; i < nRect; ++i) {
1689 // Centroid of the rectangular element.
1690 const double xc = (i + 0.5) * rectSizeX;
1691 const double yc = 0.;
1692 std::array<double, 3> centre = LocalToGlobal(xc, yc, zc, dcos, corner);
1693 // Assign element values.
1694 Element rect;
1695 rect.origin = centre;
1696 rect.lx = rectSizeX;
1697 rect.lz = elementSizeZ;
1698 rect.dA = rectSizeX * elementSizeZ;
1699 rect.dcos = dcos;
1700 // Boundary condition is applied at the origin of the rectangular
1701 // element coordinate system.
1702 rect.collocationPoint = centre;
1703 // Calculate the element vertices.
1704 const double hx = 0.5 * rectSizeX;
1705 const double hz = 0.5 * elementSizeZ;
1706 std::array<double, 3> p0 = LocalToGlobal(-hx, 0., -hz, dcos, centre);
1707 std::array<double, 3> p1 = LocalToGlobal( hx, 0., -hz, dcos, centre);
1708 std::array<double, 3> p2 = LocalToGlobal( hx, 0., hz, dcos, centre);
1709 std::array<double, 3> p3 = LocalToGlobal(-hx, 0., hz, dcos, centre);
1710 // Assign the vertices of the element.
1711 rect.xv = {p0[0], p1[0], p2[0], p3[0]};
1712 rect.yv = {p0[1], p1[1], p2[1], p3[1]};
1713 rect.zv = {p0[2], p1[2], p2[2], p3[2]};
1714 elements.push_back(std::move(rect));
1715 }
1716 }
1717 return true;
1718}
1719
1720bool ComponentNeBem3d::DiscretizeRectangle(const Primitive& primitive,
1721 const double targetSize, std::vector<Element>& elements) const {
1722
1723 std::array<double, 3> origin = {
1724 0.25 * std::accumulate(primitive.xv.begin(), primitive.xv.end(), 0.),
1725 0.25 * std::accumulate(primitive.yv.begin(), primitive.yv.end(), 0.),
1726 0.25 * std::accumulate(primitive.zv.begin(), primitive.zv.end(), 0.)};
1727
1728 // Determine the direction cosines.
1729 const double dx1 = primitive.xv[1] - primitive.xv[0];
1730 const double dy1 = primitive.yv[1] - primitive.yv[0];
1731 const double dz1 = primitive.zv[1] - primitive.zv[0];
1732 const double dx2 = primitive.xv[2] - primitive.xv[1];
1733 const double dy2 = primitive.yv[2] - primitive.yv[1];
1734 const double dz2 = primitive.zv[2] - primitive.zv[1];
1735 double lx = sqrt(dx1 * dx1 + dy1 * dy1 + dz1 * dz1);
1736 double lz = sqrt(dx2 * dx2 + dy2 * dy2 + dz2 * dz2);
1737 const std::array<double, 3> xu = {dx1 / lx, dy1 / lx, dz1 / lx};
1738 const std::array<double, 3> yu = {primitive.a, primitive.b, primitive.c};
1739 const std::array<double, 3> zu = CrossProduct(xu, yu);
1740 const std::array<std::array<double, 3>, 3> dcos = {xu, yu, zu};
1741
1742 // Determine the number of elements.
1743 unsigned int nx = NbOfSegments(lx, targetSize);
1744 unsigned int nz = NbOfSegments(lz, targetSize);
1745
1746 double elementSizeX = lx / nx;
1747 double elementSizeZ = lz / nz;
1748
1749 // Analyze the element aspect ratio.
1750 double ar = elementSizeX / elementSizeZ;
1751 if (ar > 10.) {
1752 // Reduce nz.
1753 elementSizeZ = 0.1 * elementSizeX;
1754 nz = std::max(static_cast<int>(lz / elementSizeZ), 1);
1755 elementSizeZ = lz / nz;
1756 }
1757 if (ar < 0.1) {
1758 // Reduce nx.
1759 elementSizeX = 0.1 * elementSizeZ;
1760 nx = std::max(static_cast<int>(lx / elementSizeX), 1);
1761 elementSizeX = lx / nx;
1762 }
1763
1764 const double dA = elementSizeX * elementSizeZ;
1765 for (unsigned int i = 0; i < nx; ++i) {
1766 // Centroid of the element.
1767 const double xav = -0.5 * lx + (i + 0.5) * elementSizeX;
1768 for (unsigned int k = 0; k < nz; ++k) {
1769 const double zav = -0.5 * lz + (k + 0.5) * elementSizeZ;
1770
1771 std::array<double, 3> centre = LocalToGlobal(xav, 0., zav, dcos, origin);
1772 Element element;
1773 element.origin = centre;
1774 element.lx = elementSizeX;
1775 element.lz = elementSizeZ;
1776 element.dA = dA;
1777 element.dcos = dcos;
1778 element.collocationPoint = centre;
1779 // Calculate the element vertices.
1780 const double hx = 0.5 * elementSizeX;
1781 const double hz = 0.5 * elementSizeZ;
1782 std::array<double, 3> p0 = LocalToGlobal(-hx, 0., -hz, dcos, centre);
1783 std::array<double, 3> p1 = LocalToGlobal( hx, 0., -hz, dcos, centre);
1784 std::array<double, 3> p2 = LocalToGlobal( hx, 0., hz, dcos, centre);
1785 std::array<double, 3> p3 = LocalToGlobal(-hx, 0., hz, dcos, centre);
1786 // Assign the vertices of the element.
1787 element.xv = {p0[0], p1[0], p2[0], p3[0]};
1788 element.yv = {p0[1], p1[1], p2[1], p3[1]};
1789 element.zv = {p0[2], p1[2], p2[2], p3[2]};
1790 elements.push_back(std::move(element));
1791 }
1792 }
1793 return true;
1794}
1795
1796unsigned int ComponentNeBem3d::NbOfSegments(const double length,
1797 const double target) const {
1798
1799 // Check whether the length of the primitive is long enough.
1800 if (length < MinDist) return 1;
1801 unsigned int n = static_cast<unsigned int>(length / target);
1802 if (n < m_minNbElementsOnLength) {
1803 // Need to have a minimum number of elements per primitive...
1804 n = m_minNbElementsOnLength;
1805 if (length < n * MinDist) {
1806 // ...which may not be possible if the length is small.
1807 n = static_cast<unsigned int>(length / MinDist);
1808 if (n < 1) {
1809 // However, it is necessary to have at least one element!
1810 n = 1;
1811 }
1812 }
1813 }
1814 return std::min(n, m_maxNbElementsOnLength);
1815}
1816
1817bool ComponentNeBem3d::EliminateOverlaps(const Panel& panel1,
1818 const Panel& panel2,
1819 std::vector<Panel>& panelsOut,
1820 std::vector<int>& itypo) {
1821 // *-----------------------------------------------------------------------
1822 // * PLAOVL - Isolates the parts of plane 1 that are not hidden by 2.
1823 // *-----------------------------------------------------------------------
1824
1825 const auto& xp1 = panel1.xv;
1826 const auto& yp1 = panel1.yv;
1827 const auto& zp1 = panel1.zv;
1828 const auto& xp2 = panel2.xv;
1829 const auto& yp2 = panel2.yv;
1830 const auto& zp2 = panel2.zv;
1831 // If the size of either is less than 3, simply return.
1832 if (xp1.size() <= 2 || xp2.size() <= 2) {
1833 return true;
1834 }
1835 // Compute the various tolerances.
1836 const double xmin1 = *std::min_element(std::begin(xp1), std::end(xp1));
1837 const double ymin1 = *std::min_element(std::begin(yp1), std::end(yp1));
1838 const double xmax1 = *std::max_element(std::begin(xp1), std::end(xp1));
1839 const double ymax1 = *std::max_element(std::begin(yp1), std::end(yp1));
1840
1841 const double xmin2 = *std::min_element(std::begin(xp2), std::end(xp2));
1842 const double ymin2 = *std::min_element(std::begin(yp2), std::end(yp2));
1843 const double xmax2 = *std::max_element(std::begin(xp2), std::end(xp2));
1844 const double ymax2 = *std::max_element(std::begin(yp2), std::end(yp2));
1845
1846 const double xmin = std::min(xmin1, xmin2);
1847 const double ymin = std::min(ymin1, ymin2);
1848 const double xmax = std::max(xmax1, xmax2);
1849 const double ymax = std::max(ymax1, ymax2);
1850
1851 const double epsx = 1.e-6 * std::max(std::abs(xmax), std::abs(xmin));
1852 const double epsy = 1.e-6 * std::max(std::abs(ymax), std::abs(ymin));
1853
1854 const double zsum1 = std::accumulate(std::begin(zp1), std::end(zp1), 0.);
1855 const double zsum2 = std::accumulate(std::begin(zp2), std::end(zp2), 0.);
1856 const double zmean = (zsum1 + zsum2) / (zp1.size() + zp2.size());
1857
1858 std::array<std::vector<double>, 2> xl;
1859 std::array<std::vector<double>, 2> yl;
1860 std::array<std::vector<int>, 2> flags;
1861 std::array<std::vector<double>, 2> qs;
1862 // Establish the list of special points around polygon 1.
1863 AddPoints(xp1, yp1, xp2, yp2, xl[0], yl[0], flags[0], qs[0], epsx, epsy);
1864 // Establish the list of special points around polygon 2.
1865 AddPoints(xp2, yp2, xp1, yp1, xl[1], yl[1], flags[1], qs[1], epsx, epsy);
1866
1867 bool ok = true;
1868 // Look up the cross-links: from plane 1 (2) to plane 2 (1).
1869 std::array<std::vector<int>, 2> links;
1870 for (unsigned int ic = 0; ic < 2; ++ic) {
1871 const unsigned int n1 = xl[ic].size();
1872 links[ic].assign(n1, -1);
1873 const unsigned int jc = ic == 0 ? 1 : 0;
1874 const unsigned int n2 = xl[jc].size();
1875 for (unsigned int i = 0; i < n1; ++i) {
1876 unsigned int nFound = 0;
1877 for (unsigned int j = 0; j < n2; ++j) {
1878 if (fabs(xl[ic][i] - xl[jc][j]) < epsx &&
1879 fabs(yl[ic][i] - yl[jc][j]) < epsy) {
1880 ++nFound;
1881 links[ic][i] = j;
1882 }
1883 }
1884 if (nFound == 0 && (flags[ic][i] == 2 || flags[ic][i] == 3)) {
1885 std::cerr << m_className << "::EliminateOverlaps: "
1886 << "Warning. Expected match not found (" << ic << "-"
1887 << jc << ").\n";
1888 links[ic][i] = -1;
1889 ok = false;
1890 } else if (nFound > 1) {
1891 std::cerr << m_className << "::EliminateOverlaps: "
1892 << "Warning. More than 1 match found (" << ic << "-"
1893 << jc << ").\n";
1894 links[ic][i] = -1;
1895 ok = false;
1896 }
1897 }
1898 }
1899
1900 // List the points for debugging.
1901 if (m_debug) {
1902 for (unsigned int j = 0; j < 2; ++j) {
1903 std::cout << " Polygon " << j << "\n "
1904 << " No Type x y Q links\n";
1905 const unsigned int n = xl[j].size();
1906 for (unsigned int i = 0; i < n; ++i) {
1907 printf(" %3d %5d %13.6f %13.6f %5.3f %3d\n", i, flags[j][i],
1908 xl[j][i], yl[j][i], qs[j][i], links[j][i]);
1909 }
1910 }
1911 }
1912 if (!ok) return false;
1913
1914 for (unsigned int ic = 0; ic < 2; ++ic) {
1915 // See whether all of 1 (2) is inside 2 (1).
1916 bool allInside = true;
1917 const unsigned int np = xl[ic].size();
1918 for (unsigned int i = 0; i < np; ++i) {
1919 if (flags[ic][i] != 1) {
1920 allInside = false;
1921 break;
1922 }
1923 bool inside = false, edge = false;
1924 if (ic == 0) {
1925 Polygon::Inside(xp2, yp2, xl[ic][i], yl[ic][i], inside, edge);
1926 } else {
1927 Polygon::Inside(xp1, yp1, xl[ic][i], yl[ic][i], inside, edge);
1928 }
1929 if (!(inside || edge)) {
1930 allInside = false;
1931 break;
1932 }
1933 }
1934 if (allInside) {
1935 // Apparently 1 (2) really is fully inside 2 (1).
1936 if (ic == 0) {
1937 if (m_debug) std::cout << " Curve 0 fully inside 1.\n";
1938 // Write out curve 1.
1939 panelsOut.push_back(panel1);
1940 } else {
1941 if (m_debug) std::cout << " Curve 1 fully inside 0.\n";
1942 // Write out curve 2.
1943 panelsOut.push_back(panel2);
1944 }
1945 panelsOut.back().zv.assign(panelsOut.back().xv.size(), zmean);
1946 itypo.push_back(3);
1947 std::vector<Panel> newPanels;
1948 if (ic == 0) {
1949 if (!TraceEnclosed(xl[0], yl[0], xl[1], yl[1], panel2, newPanels)) {
1950 return false;
1951 }
1952 } else {
1953 if (!TraceEnclosed(xl[1], yl[1], xl[0], yl[0], panel1, newPanels)) {
1954 return false;
1955 }
1956 }
1957 for (auto& panel : newPanels) {
1958 panel.zv.assign(panel.xv.size(), zmean);
1959 if (ic == 0) {
1960 itypo.push_back(2);
1961 } else {
1962 itypo.push_back(1);
1963 }
1964 }
1965 panelsOut.insert(panelsOut.end(), newPanels.begin(), newPanels.end());
1966 return true;
1967 }
1968 }
1969
1970 for (unsigned int ic = 0; ic < 2; ++ic) {
1971 std::vector<Panel> newPanels;
1972 const unsigned int n = xl[ic].size();
1973 // Identify the parts of 1 (2) that are not overlapped, first mark.
1974 std::vector<bool> mark(n, false);
1975 bool done = false;
1976 while (!done) {
1977 if (m_debug) {
1978 std::cout << " Searching for starting point on " << ic << ".\n";
1979 }
1980 done = true;
1981 // Try and find a new starting point
1982 for (unsigned int i = 0; i < n; ++i) {
1983 const unsigned int ii = NextPoint(i, n);
1984 // Skip parts already processed.
1985 if (mark[i] || mark[ii]) continue;
1986 // Skip if mid point is inside other volume.
1987 bool inside = false, edge = false;
1988 const double xm = 0.5 * (xl[ic][i] + xl[ic][ii]);
1989 const double ym = 0.5 * (yl[ic][i] + yl[ic][ii]);
1990 if (ic == 0) {
1991 Polygon::Inside(xp2, yp2, xm, ym, inside, edge);
1992 } else {
1993 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
1994 }
1995 if (inside || edge) continue;
1996 // Found one.
1997 done = false;
1998
1999 if (ic == 0) {
2000 // Trace this part of 1 outside 2.
2001 TraceNonOverlap(xp1, yp1, xl[0], yl[0], xl[1], yl[1], flags[0],
2002 flags[1], links[0], links[1], mark, i, panel1,
2003 newPanels);
2004 } else {
2005 // Trace this part of 2 outside 1.
2006 TraceNonOverlap(xp2, yp2, xl[1], yl[1], xl[0], yl[0], flags[1],
2007 flags[0], links[1], links[0], mark, i, panel2,
2008 newPanels);
2009 }
2010 break;
2011 }
2012 }
2013 for (auto& panel : newPanels) {
2014 panel.zv.assign(panel.xv.size(), zmean);
2015 itypo.push_back(ic + 1);
2016 }
2017 panelsOut.insert(panelsOut.end(), newPanels.begin(), newPanels.end());
2018 if (m_debug) {
2019 std::cout << " No further non-overlapped areas of " << ic << ".\n";
2020 }
2021 }
2022
2023 // Look for the overlapped parts.
2024 std::vector<Panel> newPanels;
2025 const unsigned int n1 = xl[0].size();
2026 std::vector<bool> mark1(n1, false);
2027 bool done = false;
2028 while (!done) {
2029 done = true;
2030 if (m_debug) {
2031 std::cout << " Searching for starting point on overlap.\n";
2032 }
2033 for (unsigned int i = 0; i < n1; ++i) {
2034 // Skip points already processed.
2035 if (mark1[i]) continue;
2036 // Skip if not an edge point on both 1 and 2 or internal in 2.
2037 int ip1 = i;
2038 int ip2 = links[0][ip1];
2039 if (ip2 < 0 || flags[0][ip1] == 1) {
2040 bool inside = false, edge = false;
2041 Polygon::Inside(xp2, yp2, xl[0][ip1], yl[0][ip1], inside, edge);
2042 if (!(inside || edge)) continue;
2043 } else if (flags[1][ip2] == 1) {
2044 continue;
2045 }
2046 // Found one.
2047 done = false;
2048 TraceOverlap(xp1, yp1, xp2, yp2, xl[0], yl[0], xl[1], yl[1], flags[0],
2049 links[0], links[1], mark1, ip1, ip2, panel1, newPanels);
2050 break;
2051 }
2052 }
2053 for (auto& panel : newPanels) {
2054 panel.zv.assign(panel.xv.size(), zmean);
2055 itypo.push_back(3);
2056 }
2057 panelsOut.insert(panelsOut.end(), newPanels.begin(), newPanels.end());
2058 // Finished
2059 if (m_debug) std::cout << " No further overlapped areas.\n";
2060 return true;
2061}
2062
2063bool ComponentNeBem3d::TraceEnclosed(const std::vector<double>& xl1,
2064 const std::vector<double>& yl1,
2065 const std::vector<double>& xl2,
2066 const std::vector<double>& yl2,
2067 const Panel& panel2,
2068 std::vector<Panel>& panelsOut) const {
2069 const int n1 = xl1.size();
2070 const int n2 = xl2.size();
2071 // Find 2 non-crossing connections: JP1-JP2 and KP1-KP2.
2072 unsigned int nFound = 0;
2073 int jp1 = 0, jp2 = 0;
2074 int kp1 = 0, kp2 = 0;
2075 for (int ip1 = 0; ip1 < n1; ++ip1) {
2076 const double x1 = xl1[ip1];
2077 const double y1 = yl1[ip1];
2078 for (int ip2 = 0; ip2 < n2; ++ip2) {
2079 if (nFound > 0 && ip2 == jp2) continue;
2080 const double x2 = xl2[ip2];
2081 const double y2 = yl2[ip2];
2082 bool cross = false;
2083 for (int k = 0; k < n1; ++k) {
2084 const int kk = NextPoint(k, n1);
2085 if (k == ip1 || kk == ip1) continue;
2086 double xc = 0., yc = 0.;
2087 cross =
2088 Crossing(x1, y1, x2, y2, xl1[k], yl1[k], xl1[kk], yl1[kk], xc, yc);
2089 if (cross) break;
2090 }
2091 if (cross) continue;
2092 if (m_debug) std::cout << " No crossing with 1.\n";
2093 for (int k = 0; k < n2; ++k) {
2094 const int kk = NextPoint(k, n2);
2095 if (k == ip2 || kk == ip2) continue;
2096 double xc = 0., yc = 0.;
2097 cross =
2098 Crossing(x1, y1, x2, y2, xl2[k], yl2[k], xl2[kk], yl2[kk], xc, yc);
2099 if (cross) break;
2100 }
2101 if (cross) continue;
2102 if (nFound == 0) {
2103 jp1 = ip1;
2104 jp2 = ip2;
2105 if (m_debug) {
2106 std::cout << " 1st junction: " << jp1 << ", " << jp2 << ".\n";
2107 }
2108 ++nFound;
2109 break;
2110 } else {
2111 kp1 = ip1;
2112 kp2 = ip2;
2113 double xc = 0., yc = 0.;
2114 cross = Crossing(x1, y1, x2, y2, xl1[jp1], yl1[jp1], xl2[jp2], yl2[jp2],
2115 xc, yc);
2116 if (!cross) {
2117 if (m_debug) {
2118 std::cout << " 2nd junction: " << kp1 << ", " << kp2 << ".\n";
2119 }
2120 ++nFound;
2121 break;
2122 }
2123 }
2124 }
2125 if (nFound > 1) break;
2126 }
2127 if (nFound < 2) {
2128 std::cerr << m_className << "::TraceEnclosed: Found no cut-out.\n";
2129 return false;
2130 }
2131
2132 // Create part 1 of area 2.
2133 std::vector<double> xpl;
2134 std::vector<double> ypl;
2135 if (m_debug) std::cout << " Creating part 1 of area 2.\n";
2136 for (int ip1 = jp1; ip1 <= kp1; ++ip1) {
2137 if (m_debug) std::cout << " Adding " << ip1 << " on 1.\n";
2138 xpl.push_back(xl1[ip1]);
2139 ypl.push_back(yl1[ip1]);
2140 }
2141 // Try one way.
2142 int imax = jp2 < kp2 ? jp2 + n2 : jp2;
2143 int dir = +1;
2144 for (int i = kp2; i <= imax; ++i) {
2145 int ip2 = i % n2;
2146 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
2147 xpl.push_back(xl2[ip2]);
2148 ypl.push_back(yl2[ip2]);
2149 }
2150 // Check for undesirable crossings.
2151 bool ok = true;
2152 for (int ip1 = 0; ip1 < n1; ++ip1) {
2153 if (ip1 == jp1 || ip1 == kp1) continue;
2154 bool inside = false, edge = false;
2155 Polygon::Inside(xpl, ypl, xl1[ip1], yl1[ip1], inside, edge);
2156 if (inside) {
2157 ok = false;
2158 break;
2159 }
2160 }
2161 if (!ok) {
2162 // Use the other way if this failed
2163 if (m_debug) std::cout << " Trying the other direction.\n";
2164 xpl.resize(kp1 - jp1 + 1);
2165 ypl.resize(kp1 - jp1 + 1);
2166 imax = jp2 < kp2 ? kp2 : kp2 + n2;
2167 dir = -1;
2168 for (int i = imax; i >= jp2; --i) {
2169 const int ip2 = i % n2;
2170 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
2171 xpl.push_back(xl2[ip2]);
2172 ypl.push_back(yl2[ip2]);
2173 }
2174 }
2175
2176 // Save this part.
2177 Panel newPanel1 = panel2;
2178 newPanel1.xv = xpl;
2179 newPanel1.yv = ypl;
2180 panelsOut.push_back(std::move(newPanel1));
2181 if (m_debug) std::cout << " Part 1 has " << xpl.size() << " nodes.\n";
2182
2183 // Create part 2 of area 2.
2184 xpl.clear();
2185 ypl.clear();
2186 if (m_debug) std::cout << " Creating part 2 of area 2.\n";
2187 imax = jp1 + n1;
2188 for (int i = kp1; i <= imax; ++i) {
2189 const int ip1 = i % n1;
2190 if (m_debug) std::cout << " Adding " << ip1 << " on 1.\n";
2191 xpl.push_back(xl1[ip1]);
2192 ypl.push_back(yl1[ip1]);
2193 }
2194 // Add the part over area 2.
2195 if (dir == -1) {
2196 imax = jp2 > kp2 ? jp2 : jp2 + n2;
2197 for (int i = imax; i >= kp2; --i) {
2198 const int ip2 = i % n2;
2199 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
2200 xpl.push_back(xl2[ip2]);
2201 ypl.push_back(yl2[ip2]);
2202 }
2203 } else {
2204 imax = jp2 > kp2 ? kp2 + n2 : kp2;
2205 for (int i = jp2; i <= imax; ++i) {
2206 const int ip2 = i % n2;
2207 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
2208 xpl.push_back(xl2[ip2]);
2209 ypl.push_back(yl2[ip2]);
2210 }
2211 }
2212 // Save this part.
2213 Panel newPanel2 = panel2;
2214 newPanel2.xv = xpl;
2215 newPanel2.yv = ypl;
2216 panelsOut.push_back(std::move(newPanel2));
2217 if (m_debug) std::cout << " Part 1 has " << xpl.size() << " nodes.\n";
2218 return true;
2219}
2220
2221void ComponentNeBem3d::TraceNonOverlap(
2222 const std::vector<double>& xp1, const std::vector<double>& yp1,
2223 const std::vector<double>& xl1, const std::vector<double>& yl1,
2224 const std::vector<double>& xl2, const std::vector<double>& yl2,
2225 const std::vector<int>& flags1, const std::vector<int>& flags2,
2226 const std::vector<int>& links1, const std::vector<int>& links2,
2227 std::vector<bool>& mark1, int ip1, const Panel& panel1,
2228 std::vector<Panel>& panelsOut) const {
2229 const unsigned int n1 = xl1.size();
2230 const unsigned int n2 = xl2.size();
2231
2232 // Remember the starting point.
2233 const int is1 = ip1;
2234 std::vector<double> xpl;
2235 std::vector<double> ypl;
2236 // Add the starting point.
2237 xpl.push_back(xl1[ip1]);
2238 ypl.push_back(yl1[ip1]);
2239 mark1[ip1] = true;
2240 if (m_debug) {
2241 std::cout << " Start from point " << ip1 << " on curve 1.\n";
2242 }
2243 // Next point.
2244 ip1 = NextPoint(ip1, n1);
2245 xpl.push_back(xl1[ip1]);
2246 ypl.push_back(yl1[ip1]);
2247 mark1[ip1] = true;
2248 if (m_debug) {
2249 std::cout << " Next point is " << ip1 << " on curve 1.\n";
2250 }
2251
2252 // Keep track of the curve we are currently following.
2253 unsigned int il = 1;
2254 // Direction flag (-1: backward, 0: not set, 1: forward).
2255 int dir = 0;
2256 // End-of-curve flag.
2257 bool eoc = false;
2258 int ip2 = 0;
2259 while (!eoc) {
2260 if (il == 1 && flags1[std::max(ip1, 0)] == 1) {
2261 // On curve 1 and not on the edge of curve 2?
2262 ip1 = NextPoint(ip1, n1);
2263 if (ip1 == is1) {
2264 eoc = true;
2265 continue;
2266 }
2267 mark1[ip1] = true;
2268 xpl.push_back(xl1[ip1]);
2269 ypl.push_back(yl1[ip1]);
2270 if (m_debug) {
2271 std::cout << " Went to point " << ip1 << " on curve 1.\n";
2272 }
2273 } else if (il == 1) {
2274 // On curve 1 and on the edge of curve 2?
2275 ip2 = links1[ip1];
2276 bool added = false;
2277 if (dir == +1 || dir == 0) {
2278 const double xm = 0.5 * (xl2[ip2] + xl2[NextPoint(ip2, n2)]);
2279 const double ym = 0.5 * (yl2[ip2] + yl2[NextPoint(ip2, n2)]);
2280 bool inside = false, edge = false;
2281 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2282 if (inside) {
2283 ip2 = NextPoint(ip2, n2);
2284 il = 2;
2285 dir = +1;
2286 ip1 = links2[ip2];
2287 if (ip1 == is1) {
2288 eoc = true;
2289 continue;
2290 } else if (ip1 >= 0) {
2291 mark1[ip1] = true;
2292 }
2293 xpl.push_back(xl2[ip2]);
2294 ypl.push_back(yl2[ip2]);
2295 added = true;
2296 if (m_debug) {
2297 std::cout << " Added point " << ip2 << " along 2 +.\n";
2298 }
2299 }
2300 }
2301 if (dir == -1 || dir == 0) {
2302 const double xm = 0.5 * (xl2[ip2] + xl2[PrevPoint(ip2, n2)]);
2303 const double ym = 0.5 * (yl2[ip2] + yl2[PrevPoint(ip2, n2)]);
2304 bool inside = false, edge = false;
2305 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2306 if (inside) {
2307 ip2 = PrevPoint(ip2, n2);
2308 il = 2;
2309 dir = -1;
2310 ip1 = links2[ip2];
2311 if (ip1 == is1) {
2312 eoc = true;
2313 continue;
2314 } else if (ip1 >= 0) {
2315 mark1[ip1] = true;
2316 }
2317 xpl.push_back(xl2[ip2]);
2318 ypl.push_back(yl2[ip2]);
2319 added = true;
2320 if (m_debug) {
2321 std::cout << " Added point " << ip2 << " along 2 -\n";
2322 }
2323 }
2324 }
2325 if (!added) {
2326 ip1 = NextPoint(ip1, n1);
2327 if (ip1 == is1) {
2328 eoc = true;
2329 continue;
2330 } else if (ip1 >= 0) {
2331 mark1[ip1] = true;
2332 }
2333 xpl.push_back(xl1[ip1]);
2334 ypl.push_back(yl1[ip1]);
2335 if (m_debug) std::cout << " Continued over 1.\n";
2336 }
2337 } else if (il == 2 && flags2[std::max(ip2, 0)] == 1) {
2338 // On curve 2 normal vertex (outside 1 hopefully).
2339 ip2 = dir > 0 ? NextPoint(ip2, n2) : PrevPoint(ip2, n2);
2340 ip1 = links2[ip2];
2341 if (ip1 == is1) {
2342 eoc = true;
2343 continue;
2344 } else if (ip1 >= 0) {
2345 mark1[ip1] = true;
2346 }
2347 xpl.push_back(xl2[ip2]);
2348 ypl.push_back(yl2[ip2]);
2349 if (m_debug) {
2350 std::cout << " Went to point " << ip2 << " on 2.\n";
2351 }
2352 } else if (il == 2) {
2353 // On curve 2 and on edge of 1.
2354 ip1 = links2[ip2];
2355 ip1 = NextPoint(ip1, n1);
2356 il = 1;
2357 if (ip1 == is1) {
2358 eoc = true;
2359 continue;
2360 }
2361 xpl.push_back(xl1[ip1]);
2362 ypl.push_back(yl1[ip1]);
2363 if (m_debug) std::cout << " Resumed 1 at point " << ip1 << ".\n";
2364 } else {
2365 // Other cases should not occur.
2366 std::cerr << m_className << "::TraceNonOverlap: Unexpected case.\n";
2367 return;
2368 }
2369 }
2370
2371 Panel newPanel = panel1;
2372 newPanel.xv = xpl;
2373 newPanel.yv = ypl;
2374 panelsOut.push_back(std::move(newPanel));
2375 if (m_debug) {
2376 std::cout << " End of curve reached, " << xpl.size() << " points.\n";
2377 }
2378}
2379
2380void ComponentNeBem3d::TraceOverlap(
2381 const std::vector<double>& xp1, const std::vector<double>& yp1,
2382 const std::vector<double>& xp2, const std::vector<double>& yp2,
2383 const std::vector<double>& xl1, const std::vector<double>& yl1,
2384 const std::vector<double>& xl2, const std::vector<double>& yl2,
2385 const std::vector<int>& flags1, const std::vector<int>& links1,
2386 const std::vector<int>& links2, std::vector<bool>& mark1, int ip1, int ip2,
2387 const Panel& panel1, std::vector<Panel>& panelsOut) const {
2388 int ip1L = -1;
2389 int ip1LL = -1;
2390
2391 const unsigned int n1 = xl1.size();
2392 const unsigned int n2 = xl2.size();
2393
2394 // Remember the starting points.
2395 const int is1 = ip1;
2396 const int is2 = ip2;
2397 std::vector<double> xpl;
2398 std::vector<double> ypl;
2399 xpl.push_back(xl1[ip1]);
2400 ypl.push_back(yl1[ip1]);
2401 mark1[ip1] = true;
2402
2403 // Keep track of the curve we are currently following.
2404 unsigned int il = 1;
2405 // Direction flag (-1: backward, 0: not set, 1: forward).
2406 int dir = 0;
2407 // End-of-curve flag.
2408 bool eoc = false;
2409
2410 if (m_debug) {
2411 std::cout << " Start from point " << ip1 << " on curve " << il << "\n";
2412 }
2413 while (!eoc) {
2414 ip1LL = ip1L;
2415 ip1L = ip1;
2416 if (il == 1) {
2417 // We are on curve 1. Move to the next point.
2418 const int ii = NextPoint(ip1, n1);
2419 // Maybe finished over line 1?
2420 if (ii == is1) {
2421 eoc = true;
2422 continue;
2423 }
2424 // See whether the next point of 1 is on the edge or inside of 2.
2425 bool inside = false, edge = false;
2426 if (links1[ii] >= 0) {
2427 edge = true;
2428 } else if (flags1[ii] == 1) {
2429 Polygon::Inside(xp2, yp2, xl1[ii], yl1[ii], inside, edge);
2430 }
2431 // If it is, check that it doesn't leave 2 at any stage.
2432 if (inside || edge) {
2433 const double xm = 0.5 * (xl1[ip1] + xl1[ii]);
2434 const double ym = 0.5 * (yl1[ip1] + yl1[ii]);
2435 Polygon::Inside(xp2, yp2, xm, ym, inside, edge);
2436 }
2437 // If it is, continue over 1.
2438 if (inside || edge) {
2439 ip1 = ii;
2440 if (m_debug) {
2441 std::cout << " Continued to point " << ip1 << " on "
2442 << il << "\n";
2443 }
2444 xpl.push_back(xl1[ip1]);
2445 ypl.push_back(yl1[ip1]);
2446 mark1[ip1] = true;
2447 continue;
2448 }
2449 // Else we have to continue over 2, ensure we really are on curve 2.
2450 ip2 = links1[ip1];
2451 if (ip2 < 0) {
2452 std::cerr << m_className << "::TraceOverlap: "
2453 << "No point 2 reference found; abandoned.\n";
2454 return;
2455 }
2456 // Impose a direction on 2 to avoid returning.
2457 if (dir == 0) {
2458 if (links2[NextPoint(ip2, n2)] == ip1LL &&
2459 links2[PrevPoint(ip2, n2)] == ip1LL) {
2460 std::cerr << m_className << "::TraceOverlap: "
2461 << "Both 2+ and 2- return on 1; not stored.\n";
2462 return;
2463 } else if (links2[NextPoint(ip2, n2)] == ip1LL) {
2464 if (m_debug) {
2465 std::cout << " 2+ is a return to previous point on 1.\n";
2466 }
2467 dir = -1;
2468 } else if (links2[PrevPoint(ip2, n2)] == ip1LL) {
2469 if (m_debug) {
2470 std::cout << " 2- is a return to previous point on 1.\n";
2471 }
2472 dir = +1;
2473 } else {
2474 if (m_debug) std::cout << " Both ways are OK.\n";
2475 }
2476 }
2477 // If not, try to continue over 2 in the + direction.
2478 if (dir == +1 || dir == 0) {
2479 ip2 = NextPoint(ip2, n2);
2480 if (ip2 == is2) {
2481 if (m_debug) std::cout << " Return to start over 2+.\n";
2482 eoc = true;
2483 continue;
2484 }
2485 Polygon::Inside(xp1, yp1, xl2[ip2], yl2[ip2], inside, edge);
2486 if (inside || edge) {
2487 if (m_debug) {
2488 std::cout << " Going to 2+ (point " << ip2 << " of 2).\n";
2489 }
2490 xpl.push_back(xl2[ip2]);
2491 ypl.push_back(yl2[ip2]);
2492 dir = +1;
2493 if (links2[ip2] >= 0) {
2494 ip1 = links2[ip2];
2495 mark1[ip1] = true;
2496 il = 1;
2497 if (m_debug) {
2498 std::cout << " This point is also on curve 1: "
2499 << ip1 << "\n";
2500 }
2501 } else {
2502 il = 2;
2503 }
2504 continue;
2505 }
2506 // Continuing in the + direction didn't work so go back a step.
2507 ip2 = PrevPoint(ip2, n2);
2508 }
2509 // Or if this still fails, try 2 in the - direction.
2510 if (dir == -1 || dir == 0) {
2511 ip2 = PrevPoint(ip2, n2);
2512 if (ip2 == is2) {
2513 if (m_debug) std::cout << " Return to start over 2-\n";
2514 eoc = true;
2515 continue;
2516 }
2517 Polygon::Inside(xp1, yp1, xl2[ip2], yl2[ip2], inside, edge);
2518 if (inside || edge) {
2519 if (m_debug) {
2520 std::cout << " Going to 2- (point " << ip2 << " of 2).\n";
2521 }
2522 xpl.push_back(xl2[ip2]);
2523 ypl.push_back(yl2[ip2]);
2524 dir = -1;
2525 if (links2[ip2] >= 0) {
2526 ip1 = links2[ip2];
2527 mark1[ip1] = true;
2528 il = 1;
2529 if (m_debug) {
2530 std::cout << " This point is also on 1: " << ip1 << ".\n";
2531 }
2532 } else {
2533 il = 2;
2534 }
2535 continue;
2536 }
2537 }
2538 // Should not get here.
2539 if (m_debug) std::cout << " Dead end.\n";
2540 return;
2541 } else if (il == 2) {
2542 // We are on curve 2. Ensure the direction is set.
2543 if (dir == 0) {
2544 std::cerr << m_className << "::TraceOverlap: "
2545 << "Direction not set; abandoned.\n";
2546 return;
2547 }
2548 // Move to the next point.
2549 ip2 = dir > 0 ? NextPoint(ip2, n2) : PrevPoint(ip2, n2);
2550 // Maybe finished over line 2?
2551 if (ip2 == is2) {
2552 // Reached the end.
2553 eoc = true;
2554 continue;
2555 }
2556 // Next step over 2.
2557 if (m_debug) {
2558 std::cout << " Stepped over 2 to point " << ip2 << " of 2.\n";
2559 }
2560 xpl.push_back(xl2[ip2]);
2561 ypl.push_back(yl2[ip2]);
2562 if (links2[ip2] >= 0) {
2563 ip1 = links2[ip2];
2564 mark1[ip1] = true;
2565 il = 1;
2566 if (m_debug) {
2567 std::cout << " This point is also on curve 1: " << ip1 << ".\n";
2568 }
2569 } else {
2570 il = 2;
2571 }
2572 }
2573 }
2574
2575 if (xpl.size() <= 2) {
2576 if (m_debug) std::cout << " Too few points.\n";
2577 } else {
2578 Panel newPanel = panel1;
2579 newPanel.xv = xpl;
2580 newPanel.yv = ypl;
2581 panelsOut.push_back(std::move(newPanel));
2582 }
2583
2584 if (m_debug) {
2585 std::cout << " End of curve reached, " << xpl.size() << " points.\n";
2586 }
2587}
2588
2589bool ComponentNeBem3d::MakePrimitives(const Panel& panelIn,
2590 std::vector<Panel>& panelsOut) const {
2591 // *-----------------------------------------------------------------------
2592 // * PLATRC - Cuts a polygon into right-angled triangles.
2593 // *-----------------------------------------------------------------------
2594
2595 // Establish tolerances.
2596 // TODO! Class member?
2597 const double epsang = 1.e-6; // BEMEPA
2598
2599 if (panelIn.xv.empty() || panelIn.yv.empty() || panelIn.zv.empty()) {
2600 return false;
2601 }
2602 // Determine the mean z value.
2603 const double zsum =
2604 std::accumulate(std::begin(panelIn.zv), std::end(panelIn.zv), 0.);
2605 const double zmean = zsum / panelIn.zv.size();
2606
2607 std::vector<Panel> stack;
2608 stack.push_back(panelIn);
2609 stack.back().zv.clear();
2610 for (unsigned int k = 0; k < stack.size(); ++k) {
2611 // Next polygon.
2612 const auto& xp1 = stack[k].xv;
2613 const auto& yp1 = stack[k].yv;
2614 const unsigned int np = xp1.size();
2615 if (m_debug) {
2616 std::cout << " Polygon " << k << " with " << np << " nodes.\n";
2617 }
2618 if (np <= 2) {
2619 // Too few nodes.
2620 if (m_debug) std::cout << " Too few points.\n";
2621 continue;
2622 }
2623
2624 // See whether this is a right-angled triangle.
2625 if (np == 3) {
2626 const double x12 = xp1[0] - xp1[1];
2627 const double y12 = yp1[0] - yp1[1];
2628 const double x23 = xp1[1] - xp1[2];
2629 const double y23 = yp1[1] - yp1[2];
2630 const double x31 = xp1[2] - xp1[0];
2631 const double y31 = yp1[2] - yp1[0];
2632 const double x32 = xp1[2] - xp1[1];
2633 const double y32 = yp1[2] - yp1[1];
2634 const double x13 = xp1[0] - xp1[2];
2635 const double y13 = yp1[0] - yp1[2];
2636 const double x21 = xp1[1] - xp1[0];
2637 const double y21 = yp1[1] - yp1[0];
2638 const double s12 = x12 * x12 + y12 * y12;
2639 const double s32 = x32 * x32 + y32 * y32;
2640 const double s13 = x13 * x13 + y13 * y13;
2641 const double s23 = x23 * x23 + y23 * y23;
2642 const double s31 = x31 * x31 + y31 * y31;
2643 const double s21 = x21 * x21 + y21 * y21;
2644 if (fabs(x12 * x32 + y12 * y32) < epsang * sqrt(s12 * s32)) {
2645 if (m_debug) {
2646 std::cout << " Right-angled triangle node 2 - done.\n";
2647 }
2648 panelsOut.push_back(stack[k]);
2649 continue;
2650 } else if (fabs(x13 * x23 + y13 * y23) < epsang * sqrt(s13 * s23)) {
2651 if (m_debug) {
2652 std::cout << " Right-angled triangle node 3 - rearrange.\n";
2653 }
2654 Panel panel = stack[k];
2655 panel.xv = {xp1[1], xp1[2], xp1[0]};
2656 panel.yv = {yp1[1], yp1[2], yp1[0]};
2657 panelsOut.push_back(std::move(panel));
2658 continue;
2659 } else if (fabs(x31 * x21 + y31 * y21) < epsang * sqrt(s31 * s21)) {
2660 if (m_debug) {
2661 std::cout << " Right-angled triangle node 1 - rearrange.\n";
2662 }
2663 Panel panel = stack[k];
2664 panel.xv = {xp1[2], xp1[0], xp1[1]};
2665 panel.yv = {yp1[2], yp1[0], yp1[1]};
2666 panelsOut.push_back(std::move(panel));
2667 continue;
2668 }
2669 }
2670 // See whether this is a rectangle.
2671 if (np == 4) {
2672 const double x12 = xp1[0] - xp1[1];
2673 const double y12 = yp1[0] - yp1[1];
2674 const double x23 = xp1[1] - xp1[2];
2675 const double y23 = yp1[1] - yp1[2];
2676 const double x34 = xp1[2] - xp1[3];
2677 const double y34 = yp1[2] - yp1[3];
2678 const double x43 = xp1[3] - xp1[2];
2679 const double y43 = yp1[3] - yp1[2];
2680 const double x14 = xp1[0] - xp1[3];
2681 const double y14 = yp1[0] - yp1[3];
2682 const double x32 = xp1[2] - xp1[1];
2683 const double y32 = yp1[2] - yp1[1];
2684
2685 const double s12 = x12 * x12 + y12 * y12;
2686 const double s23 = x23 * x23 + y23 * y23;
2687 const double s32 = x32 * x32 + y32 * y32;
2688 const double s34 = x34 * x34 + y34 * y34;
2689 const double s43 = x43 * x43 + y43 * y43;
2690 const double s14 = x14 * x14 + y14 * y14;
2691 if (fabs(x12 * x32 + y12 * y32) < epsang * sqrt(s12 * s32) &&
2692 fabs(x23 * x43 + y23 * y43) < epsang * sqrt(s23 * s43) &&
2693 fabs(x14 * x34 + y14 * y34) < epsang * sqrt(s14 * s34)) {
2694 if (m_debug) std::cout << " Rectangle.\n";
2695 panelsOut.push_back(stack[k]);
2696 continue;
2697 }
2698 }
2699 // See whether there are parallel sides, e.g. a trapezium (UK English).
2700 if (np >= 4 && SplitTrapezium(stack[k], stack, panelsOut, epsang)) continue;
2701
2702 // Find a right-angled corner we can cut off.
2703 if (m_debug) std::cout << " Trying to find a right-angle\n";
2704 bool corner = false;
2705 for (unsigned int ip = 0; ip < np; ++ip) {
2706 // Take only right angles.
2707 const unsigned int inext = NextPoint(ip, np);
2708 const unsigned int iprev = PrevPoint(ip, np);
2709
2710 const double dxprev = xp1[iprev] - xp1[ip];
2711 const double dyprev = yp1[iprev] - yp1[ip];
2712 const double dxnext = xp1[inext] - xp1[ip];
2713 const double dynext = yp1[inext] - yp1[ip];
2714 if (fabs(dxprev * dxnext + dyprev * dynext) >
2715 epsang * sqrt((dxprev * dxprev + dyprev * dyprev) *
2716 (dxnext * dxnext + dynext * dynext))) {
2717 continue;
2718 }
2719 // Ensure the midpoint is internal.
2720 if (np > 3) {
2721 const double xm = 0.5 * (xp1[iprev] + xp1[inext]);
2722 const double ym = 0.5 * (yp1[iprev] + yp1[inext]);
2723 bool inside = false, edge = false;
2724 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2725 if (!inside) continue;
2726 }
2727 // Check all vertex crossings.
2728 bool cross = false;
2729 for (unsigned int jp = 0; jp < np; ++jp) {
2730 // Accept immediate contact.
2731 const unsigned int jnext = NextPoint(jp, np);
2732 if (jp == iprev || jp == ip || jp == inext || jnext == iprev ||
2733 jnext == ip || jnext == inext)
2734 continue;
2735 // Check crossing.
2736 double xc = 0., yc = 0.;
2737 if (Crossing(xp1[iprev], yp1[iprev], xp1[inext], yp1[inext], xp1[jp],
2738 yp1[jp], xp1[jnext], yp1[jnext], xc, yc)) {
2739 cross = true;
2740 break;
2741 }
2742 }
2743 if (cross) continue;
2744 // Found a triangle, introduce shorthand node references.
2745 if (m_debug) {
2746 std::cout << " Cutting at right-angled corner " << ip << "\n";
2747 }
2748 corner = true;
2749 Panel panel = stack[k];
2750 panel.xv = {xp1[iprev], xp1[ip], xp1[inext]};
2751 panel.yv = {yp1[iprev], yp1[ip], yp1[inext]};
2752 stack.push_back(std::move(panel));
2753 // Eliminate this node from the polygon.
2754 stack[k].xv.erase(stack[k].xv.begin() + ip);
2755 stack[k].yv.erase(stack[k].yv.begin() + ip);
2756 stack.push_back(std::move(stack[k]));
2757 if (m_debug) {
2758 std::cout << " Going for another pass, NP = " << np << "\n";
2759 }
2760 break;
2761 }
2762 if (corner) continue;
2763
2764 // Find any corner we can cut off.
2765 if (m_debug) std::cout << " Trying to find a corner\n";
2766 corner = false;
2767 for (unsigned int ip = 0; ip < np; ++ip) { // 20
2768 const unsigned int iprev = PrevPoint(ip, np);
2769 const unsigned int inext = NextPoint(ip, np);
2770 // Ensure the midpoint is internal.
2771 if (np > 3) {
2772 const double xm = 0.5 * (xp1[iprev] + xp1[inext]);
2773 const double ym = 0.5 * (yp1[iprev] + yp1[inext]);
2774 bool inside = false, edge = false;
2775 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2776 if (!inside) continue;
2777 }
2778 // Check all vertex crossings.
2779 bool cross = false;
2780 for (unsigned int jp = 0; jp < np; ++jp) {
2781 const unsigned int jj = NextPoint(jp, np);
2782 // Accept immediate contact.
2783 if (jp == iprev || jp == ip || jp == inext || jj == iprev || jj == ip ||
2784 jj == inext)
2785 continue;
2786 // Check crossing.
2787 double xc = 0., yc = 0.;
2788 if (Crossing(xp1[iprev], yp1[iprev], xp1[inext], yp1[inext], xp1[jp],
2789 yp1[jp], xp1[jj], yp1[jj], xc, yc)) {
2790 cross = true;
2791 break;
2792 }
2793 }
2794 if (cross) continue;
2795 // Found a triangle, introduce shorthand node references.
2796 if (m_debug) std::cout << " Cutting at corner " << ip << "\n";
2797 corner = true;
2798 const double x1 = xp1[iprev];
2799 const double x2 = xp1[ip];
2800 const double x3 = xp1[inext];
2801 const double y1 = yp1[iprev];
2802 const double y2 = yp1[ip];
2803 const double y3 = yp1[inext];
2804 const double s21 = (x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1);
2805 const double s31 = (x3 - x1) * (x3 - x1) + (y3 - y1) * (y3 - y1);
2806 const double s32 = (x3 - x2) * (x3 - x2) + (y3 - y2) * (y3 - y2);
2807 const double s12 = (x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2);
2808 const double s13 = (x1 - x3) * (x1 - x3) + (y1 - y3) * (y1 - y3);
2809 const double s23 = (x2 - x3) * (x2 - x3) + (y2 - y3) * (y2 - y3);
2810 // Find the biggest opening angle.
2811 const double a1 =
2812 ((x2 - x1) * (x3 - x1) + (y2 - y1) * (y3 - y1)) / sqrt(s21 * s31);
2813 const double a2 =
2814 ((x3 - x2) * (x1 - x2) + (y3 - y2) * (y1 - y2)) / sqrt(s32 * s12);
2815 const double a3 =
2816 ((x1 - x3) * (x2 - x3) + (y1 - y3) * (y2 - y3)) / sqrt(s13 * s23);
2817 if (m_debug) {
2818 const double phi1 = acos(a1);
2819 const double phi2 = acos(a2);
2820 const double phi3 = acos(a3);
2821 std::cout << " Angles: " << RadToDegree * phi1 << ", "
2822 << RadToDegree * phi2 << ", " << RadToDegree * phi3 << "\n";
2823 const double sumphi = phi1 + phi2 + phi3;
2824 std::cout << " Sum = " << RadToDegree * sumphi << "\n";
2825 }
2826 // See whether one angle is more or less right-angled.
2827 if (fabs(a1) < epsang || fabs(a2) < epsang || fabs(a3) < epsang) {
2828 if (m_debug) std::cout << " Right-angled corner cut off.\n";
2829 Panel panel = stack[k];
2830 if (fabs(a1) < epsang) {
2831 panel.xv = {x3, x1, x2};
2832 panel.yv = {y3, y1, y2};
2833 } else if (fabs(a2) < epsang) {
2834 panel.xv = {x1, x2, x3};
2835 panel.yv = {y1, y2, y3};
2836 } else {
2837 panel.xv = {x2, x3, x1};
2838 panel.yv = {y2, y3, y1};
2839 }
2840 stack.push_back(std::move(panel));
2841 } else if (a1 <= a2 && a1 <= a3) {
2842 if (m_debug) std::cout << " A1 < A2, A3 - adding 2 triangles.\n";
2843 const double xc = x2 + a2 * (x3 - x2) * sqrt(s12 / s32);
2844 const double yc = y2 + a2 * (y3 - y2) * sqrt(s12 / s32);
2845 Panel panel1 = stack[k];
2846 panel1.xv = {x3, xc, x1};
2847 panel1.yv = {y3, yc, y1};
2848 stack.push_back(std::move(panel1));
2849 Panel panel2 = stack[k];
2850 panel2.xv = {x2, xc, x1};
2851 panel2.yv = {y2, yc, y1};
2852 stack.push_back(std::move(panel2));
2853 } else if (a2 <= a1 && a2 <= a3) {
2854 if (m_debug) std::cout << " A2 < A1, A3 - adding 2 triangles.\n";
2855 const double xc = x3 + a3 * (x1 - x3) * sqrt(s23 / s13);
2856 const double yc = y3 + a3 * (y1 - y3) * sqrt(s23 / s13);
2857 Panel panel1 = stack[k];
2858 panel1.xv = {x1, xc, x2};
2859 panel1.yv = {y1, yc, y2};
2860 stack.push_back(std::move(panel1));
2861 Panel panel2 = stack[k];
2862 panel2.xv = {x3, xc, x2};
2863 panel2.yv = {y3, yc, y2};
2864 stack.push_back(std::move(panel2));
2865 } else {
2866 if (m_debug) std::cout << " A3 < A1, A2 - adding 2 triangles.\n";
2867 const double xc = x1 + a1 * (x2 - x1) * sqrt(s31 / s21);
2868 const double yc = y1 + a1 * (y2 - y1) * sqrt(s31 / s21);
2869 Panel panel1 = stack[k];
2870 panel1.xv = {x1, xc, x3};
2871 panel1.yv = {y1, yc, y3};
2872 stack.push_back(std::move(panel1));
2873 Panel panel2 = stack[k];
2874 panel2.xv = {x2, xc, x3};
2875 panel2.yv = {y2, yc, y3};
2876 stack.push_back(std::move(panel2));
2877 }
2878 // Eliminate this node from the polygon.
2879 stack[k].xv.erase(stack[k].xv.begin() + ip);
2880 stack[k].yv.erase(stack[k].yv.begin() + ip);
2881 stack.push_back(std::move(stack[k]));
2882 if (m_debug) {
2883 std::cout << " Going for another pass, NP = " << np << "\n";
2884 }
2885 break;
2886 }
2887 if (corner) continue;
2888 std::cerr << m_className << "::MakePrimitives:\n "
2889 << "Unable to identify a corner to cut, probably a degenerate "
2890 "polygon.\n";
2891 // Next stack element.
2892 }
2893 for (auto& panel : panelsOut) {
2894 panel.zv.assign(panel.xv.size(), zmean);
2895 }
2896 return true;
2897}
2898
2899bool ComponentNeBem3d::SplitTrapezium(const Panel panelIn,
2900 std::vector<Panel>& stack,
2901 std::vector<Panel>& panelsOut,
2902 const double epsang) const {
2903 const auto xp1 = panelIn.xv;
2904 const auto yp1 = panelIn.yv;
2905 const unsigned int np = xp1.size();
2906 for (unsigned int ip = 0; ip < np; ++ip) {
2907 const unsigned int inext = NextPoint(ip, np);
2908 const double xi0 = xp1[ip];
2909 const double yi0 = yp1[ip];
2910 const double xi1 = xp1[inext];
2911 const double yi1 = yp1[inext];
2912 const double dxi = xi0 - xi1;
2913 const double dyi = yi0 - yi1;
2914 const double si2 = dxi * dxi + dyi * dyi;
2915 for (unsigned int jp = ip + 2; jp < np; ++jp) {
2916 const unsigned int jnext = NextPoint(jp, np);
2917 // Skip adjacent segments.
2918 if (ip == jp || ip == jnext || inext == jp || inext == jnext) {
2919 continue;
2920 }
2921 const double xj0 = xp1[jp];
2922 const double yj0 = yp1[jp];
2923 const double xj1 = xp1[jnext];
2924 const double yj1 = yp1[jnext];
2925 const double dxj = xj0 - xj1;
2926 const double dyj = yj0 - yj1;
2927 const double sj2 = dxj * dxj + dyj * dyj;
2928 // Require parallelism.
2929 const double sij = sqrt(si2 * sj2);
2930 if (fabs(dxi * dxj + dyi * dyj + sij) > epsang * sij) continue;
2931 if (m_debug) {
2932 std::cout << " Found parallel sections: "
2933 << ip << ", " << jp << "\n";
2934 }
2935 // Avoid division by zero.
2936 if (sj2 <= 0 || si2 <= 0) {
2937 std::cerr << m_className << "::SplitTrapezium:\n "
2938 << "Zero norm segment found; skipped.\n";
2939 continue;
2940 }
2941 // Establish the cutting lines.
2942 const double xl1 =
2943 ((xi0 - xj0) * (xj1 - xj0) + (yi0 - yj0) * (yj1 - yj0)) / sj2;
2944 const double xl2 =
2945 ((xi1 - xj0) * (xj1 - xj0) + (yi1 - yj0) * (yj1 - yj0)) / sj2;
2946 const double xl3 =
2947 ((xj0 - xi0) * (xi1 - xi0) + (yj0 - yi0) * (yi1 - yi0)) / si2;
2948 const double xl4 =
2949 ((xj1 - xi0) * (xi1 - xi0) + (yj1 - yi0) * (yi1 - yi0)) / si2;
2950 if (m_debug) {
2951 std::cout << " xl1 = " << xl1 << ", xl2 = " << xl2 << ", "
2952 << "xl3 = " << xl3 << ", xl4 = " << xl4 << "\n";
2953 }
2954 // Check that there is at all a rectangle.
2955 const double r1 = (xl1 + epsang) * (1. + epsang - xl1);
2956 const double r2 = (xl2 + epsang) * (1. + epsang - xl2);
2957 const double r3 = (xl3 + epsang) * (1. + epsang - xl3);
2958 const double r4 = (xl4 + epsang) * (1. + epsang - xl4);
2959 if ((r1 < 0 && r4 < 0) || (r2 < 0 && r3 < 0)) {
2960 if (m_debug) std::cout << " No rectangle.\n";
2961 continue;
2962 }
2963 // Determine the rectangular part.
2964 std::vector<double> xpl(4, 0.);
2965 std::vector<double> ypl(4, 0.);
2966 if (r1 >= 0) {
2967 xpl[0] = xi0;
2968 ypl[0] = yi0;
2969 xpl[1] = xj0 + xl1 * (xj1 - xj0);
2970 ypl[1] = yj0 + xl1 * (yj1 - yj0);
2971 } else if (r4 >= 0) {
2972 xpl[0] = xi0 + xl4 * (xi1 - xi0);
2973 ypl[0] = yi0 + xl4 * (yi1 - yi0);
2974 xpl[1] = xj1;
2975 ypl[1] = yj1;
2976 }
2977 if (r2 >= 0) {
2978 xpl[2] = xj0 + xl2 * (xj1 - xj0);
2979 ypl[2] = yj0 + xl2 * (yj1 - yj0);
2980 xpl[3] = xi1;
2981 ypl[3] = yi1;
2982 } else if (r3 >= 0) {
2983 xpl[2] = xj0;
2984 ypl[2] = yj0;
2985 xpl[3] = xi0 + xl3 * (xi1 - xi0);
2986 ypl[3] = yi0 + xl3 * (yi1 - yi0);
2987 }
2988 // Verify that the midpoints of these lines are internal.
2989 double xm = 0.5 * (xpl[0] + xpl[1]);
2990 double ym = 0.5 * (ypl[0] + ypl[1]);
2991 bool inside = false, edge = false;
2992 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2993 if (!(inside || edge)) {
2994 if (m_debug) std::cout << " Midpoint 1 not internal.\n";
2995 continue;
2996 }
2997 xm = 0.5 * (xpl[2] + xpl[3]);
2998 ym = 0.5 * (ypl[2] + ypl[3]);
2999 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
3000 if (!(inside || edge)) {
3001 if (m_debug) std::cout << " Midpoint 2 not internal.\n";
3002 continue;
3003 }
3004
3005 const unsigned int iprev = PrevPoint(ip, np);
3006 const unsigned int jprev = PrevPoint(jp, np);
3007 // Ensure there are no crossings, accepting contact.
3008 bool cross = false;
3009 for (unsigned int i = 0; i < np; ++i) {
3010 if ((i == iprev && r1 >= 0) || i == ip || (i == inext && r2 >= 0) ||
3011 (i == jprev && r3 >= 0) || i == jp || (i == jnext && r4 >= 0))
3012 continue;
3013 const unsigned int ii = NextPoint(i, np);
3014 double xc = 0., yc = 0.;
3015 if (Crossing(xp1[i], yp1[i], xp1[ii], yp1[ii], xpl[0], ypl[0], xpl[1],
3016 ypl[1], xc, yc) ||
3017 Crossing(xp1[i], yp1[i], xp1[ii], yp1[ii], xpl[2], ypl[2], xpl[3],
3018 ypl[3], xc, yc)) {
3019 if (m_debug) {
3020 std::cout << " Crossing (edge " << i << ", " << ii
3021 << "), ip = " << ip << ", jp = " << jp << "\n";
3022 }
3023 cross = true;
3024 break;
3025 }
3026 }
3027 if (cross) continue;
3028 // Add the rectangular part.
3029 if ((fabs(xl1) < epsang && fabs(xl3) < epsang) ||
3030 (fabs(1. - xl2) < epsang && fabs(1. - xl4) < epsang)) {
3031 if (m_debug) std::cout << " Not stored, degenerate.\n";
3032 } else {
3033 if (m_debug) std::cout << " Adding rectangle.\n";
3034 Panel panel = panelIn;
3035 panel.xv = xpl;
3036 panel.yv = ypl;
3037 panelsOut.push_back(std::move(panel));
3038 }
3039 // First non-rectangular section.
3040 xpl.clear();
3041 ypl.clear();
3042 if (m_debug) std::cout << " First non-rectangular section.\n";
3043 for (unsigned int i = jp + 1; i <= ip + np; ++i) {
3044 const unsigned int ii = i % np;
3045 xpl.push_back(xp1[ii]);
3046 ypl.push_back(yp1[ii]);
3047 }
3048 if (r1 >= 0 && r4 >= 0) {
3049 if (m_debug) std::cout << " 1-4 degenerate\n";
3050 } else if (r1 >= 0) {
3051 if (m_debug) std::cout << " Using 1\n";
3052 xpl.push_back(xj0 + xl1 * (xj1 - xj0));
3053 ypl.push_back(yj0 + xl1 * (yj1 - yj0));
3054 } else if (r4 >= 0) {
3055 if (m_debug) std::cout << " Using 4\n";
3056 xpl.push_back(xi0 + xl4 * (xi1 - xi0));
3057 ypl.push_back(yi0 + xl4 * (yi1 - yi0));
3058 } else {
3059 if (m_debug) std::cout << " Neither 1 nor 4, should not happen\n";
3060 }
3061 if (xpl.size() < 3) {
3062 if (m_debug) {
3063 std::cout << " Not stored, only " << xpl.size()
3064 << " vertices.\n";
3065 }
3066 } else {
3067 if (m_debug) {
3068 std::cout << " Adding non-rectangular part with " << xpl.size()
3069 << " vertices.\n";
3070 }
3071 Panel panel = panelIn;
3072 panel.xv = xpl;
3073 panel.yv = ypl;
3074 stack.push_back(std::move(panel));
3075 }
3076 // Second non-rectangular section.
3077 xpl.clear();
3078 ypl.clear();
3079 if (m_debug) std::cout << " Second non-rectangular section.\n";
3080 for (unsigned int i = ip + 1; i <= jp; ++i) {
3081 const unsigned int ii = i % np;
3082 xpl.push_back(xp1[ii]);
3083 ypl.push_back(yp1[ii]);
3084 }
3085 if (r2 >= 0 && r3 >= 0) {
3086 if (m_debug) std::cout << " 2-3 degenerate\n";
3087 } else if (r2 >= 0) {
3088 if (m_debug) std::cout << " Using 2\n";
3089 xpl.push_back(xj0 + xl2 * (xj1 - xj0));
3090 ypl.push_back(yj0 + xl2 * (yj1 - yj0));
3091 } else if (r3 >= 0) {
3092 if (m_debug) std::cout << " Using 3\n";
3093 xpl.push_back(xi0 + xl3 * (xi1 - xi0));
3094 ypl.push_back(yi0 + xl3 * (yi1 - yi0));
3095 } else {
3096 if (m_debug) std::cout << " Neither 2 nor 3, should not happen\n";
3097 }
3098 if (xpl.size() < 3) {
3099 if (m_debug) {
3100 std::cout << " Not stored, only " << xpl.size()
3101 << " vertices.\n";
3102 }
3103 } else {
3104 if (m_debug) {
3105 std::cout << " Adding non-rectangular part with " << xpl.size()
3106 << " vertices.\n";
3107 }
3108 Panel panel = panelIn;
3109 panel.xv = xpl;
3110 panel.yv = ypl;
3111 stack.push_back(std::move(panel));
3112 }
3113 return true;
3114 }
3115 }
3116 return false;
3117}
3118
3119bool ComponentNeBem3d::GetPrimitive(const unsigned int i, double& a, double& b,
3120 double& c, std::vector<double>& xv,
3121 std::vector<double>& yv,
3122 std::vector<double>& zv, int& interface,
3123 double& v, double& q,
3124 double& lambda) const {
3125 if (i >= m_primitives.size()) {
3126 std::cerr << m_className << "::GetPrimitive: Index out of range.\n";
3127 return false;
3128 }
3129 const auto& primitive = m_primitives[i];
3130 a = primitive.a;
3131 b = primitive.b;
3132 c = primitive.c;
3133 xv = primitive.xv;
3134 yv = primitive.yv;
3135 zv = primitive.zv;
3136 interface = primitive.interface;
3137 v = primitive.v;
3138 q = primitive.q;
3139 lambda = primitive.lambda;
3140 return true;
3141}
3142
3143bool ComponentNeBem3d::GetPrimitive(const unsigned int i, double& a, double& b,
3144 double& c, std::vector<double>& xv,
3145 std::vector<double>& yv,
3146 std::vector<double>& zv, int& vol1, int& vol2) const {
3147 if (i >= m_primitives.size()) {
3148 std::cerr << m_className << "::GetPrimitive: Index out of range.\n";
3149 return false;
3150 }
3151 const auto& primitive = m_primitives[i];
3152 a = primitive.a;
3153 b = primitive.b;
3154 c = primitive.c;
3155 xv = primitive.xv;
3156 yv = primitive.yv;
3157 zv = primitive.zv;
3158 vol1 = primitive.vol1;
3159 vol2 = primitive.vol2;
3160 return true;
3161}
3162
3163bool ComponentNeBem3d::GetVolume(const unsigned int vol, int& shape,
3164 int& material, double& epsilon,
3165 double& potential, double& charge,
3166 int& bc) {
3167
3168 if (!m_geometry) return false;
3169 const unsigned int nSolids = m_geometry->GetNumberOfSolids();
3170 for (unsigned int i = 0; i < nSolids; ++i) {
3171 Medium* medium = nullptr;
3172 const auto solid = m_geometry->GetSolid(i, medium);
3173 if (!solid) continue;
3174 if (solid->GetId() != vol) continue;
3175 if (solid->IsTube() || solid->IsWire()) {
3176 shape = 1;
3177 } else if (solid->IsHole()) {
3178 shape = 2;
3179 } else if (solid->IsBox()) {
3180 shape = 3;
3181 } else if (solid->IsSphere()) {
3182 shape = 4;
3183 } else if (solid->IsRidge()) {
3184 shape = 5;
3185 } else if (solid->IsExtrusion()) {
3186 shape = 6;
3187 } else {
3188 std::cerr << m_className << "::GetVolume: Unknown solid shape.\n";
3189 return false;
3190 }
3191 material = medium ? medium->GetId() : 11;
3192 epsilon = medium ? medium->GetDielectricConstant() : 1.;
3193 potential = solid->GetBoundaryPotential();
3194 charge = solid->GetBoundaryChargeDensity();
3195 bc = solid->GetBoundaryConditionType();
3196 return true;
3197 }
3198 return false;
3199}
3200
3201int ComponentNeBem3d::GetVolume(const double x, const double y,
3202 const double z) {
3203 if (!m_geometry) return -1;
3204 const unsigned int nSolids = m_geometry->GetNumberOfSolids();
3205 for (unsigned int i = 0; i < nSolids; ++i) {
3206 Medium* medium = nullptr;
3207 const auto solid = m_geometry->GetSolid(i, medium);
3208 if (!solid) continue;
3209 if (solid->IsInside(x, y, z)) return solid->GetId();
3210 }
3211 return -1;
3212}
3213
3214bool ComponentNeBem3d::GetElement(const unsigned int i,
3215 std::vector<double>& xv,
3216 std::vector<double>& yv,
3217 std::vector<double>& zv, int& interface,
3218 double& bc, double& lambda) const {
3219
3220 if (i >= m_elements.size()) {
3221 std::cerr << m_className << "::GetElement: Index out of range.\n";
3222 return false;
3223 }
3224 const auto& element = m_elements[i];
3225 xv = element.xv;
3226 yv = element.yv;
3227 zv = element.zv;
3228 interface = element.interface;
3229 bc = element.bc;
3230 lambda = element.lambda;
3231 return true;
3232}
3233
3235 m_primitives.clear();
3236 m_elements.clear();
3237 m_ynplan.fill(false);
3238 m_coplan.fill(0.);
3239 m_vtplan.fill(0.);
3240 m_ready = false;
3241}
3242
3244
3245 for (unsigned int i = 0; i < 3; ++i) {
3246 // Cannot have regular and mirror periodicity at the same time.
3247 if (m_periodic[i] && m_mirrorPeriodic[i]) {
3248 std::cerr << m_className << "::UpdatePeriodicity:\n"
3249 << " Both simple and mirror periodicity requested. Reset.\n";
3250 m_periodic[i] = false;
3251 m_mirrorPeriodic[i] = false;
3252 continue;
3253 }
3254 if ((m_periodic[i] || m_mirrorPeriodic[i]) && m_periodicLength[i] < Small) {
3255 std::cerr << m_className << "::UpdatePeriodicity:\n"
3256 << " Periodic length is not set. Reset.\n";
3257 m_periodic[i] = m_mirrorPeriodic[i] = false;
3258 }
3259 }
3260
3262 std::cerr << m_className << "::UpdatePeriodicity:\n"
3263 << " Axial periodicity is not available.\n";
3264 m_axiallyPeriodic.fill(false);
3265 }
3266
3269 std::cerr << m_className << "::UpdatePeriodicity:\n"
3270 << " Rotation symmetry is not available.\n";
3271 m_rotationSymmetric.fill(false);
3272 }
3273
3274}
3275}
bool GetPlaneZ(const unsigned int i, double &z, double &v) const
Retrieve the parameters of a plane at constant z.
bool GetPlaneY(const unsigned int i, double &y, double &v) const
Retrieve the parameters of a plane at constant y.
void AddPlaneZ(const double z, const double voltage)
Add a plane at constant z.
bool GetPeriodicityX(double &s) const
Get the periodic length in the x-direction.
bool GetPeriodicityZ(double &s) const
Get the periodic length in the z-direction.
void Reset() override
Reset the component.
void SetMirrorPeriodicityX(const double s)
Set the periodic length [cm] in the x-direction.
void SetPeriodicityX(const double s)
Set the periodic length [cm] in the x-direction.
void UpdatePeriodicity() override
Verify periodicities.
bool GetVoltageRange(double &vmin, double &vmax) override
Calculate the voltage range [V].
unsigned int GetNumberOfPlanesY() const
Get the number of equipotential planes at constant y.
bool GetPeriodicityY(double &s) const
Get the periodic length in the y-direction.
unsigned int GetNumberOfPlanesZ() const
Get the number of equipotential planes at constant z.
void AddPlaneY(const double y, const double voltage)
Add a plane at constant y.
void AddPlaneX(const double x, const double voltage)
Add a plane at constant x.
void SetMirrorPeriodicityZ(const double s)
Set the periodic length [cm] in the z-direction.
void SetPeriodicityZ(const double s)
Set the periodic length [cm] in the z-direction.
bool GetVolume(const unsigned int vol, int &shape, int &material, double &eps, double &potential, double &charge, int &bc)
void SetPeriodicityY(const double s)
Set the periodic length [cm] in the y-direction.
Medium * GetMedium(const double x, const double y, const double z) override
Get the medium at a given location (x, y, z).
unsigned int GetNumberOfPlanesX() const
Get the number of equipotential planes at constant x.
bool GetElement(const unsigned int i, std::vector< double > &xv, std::vector< double > &yv, std::vector< double > &zv, int &interface, double &bc, double &lambda) const
void SetTargetElementSize(const double length)
void SetMirrorPeriodicityY(const double s)
Set the periodic length [cm] in the y-direction.
void SetPeriodicCopies(const unsigned int nx, const unsigned int ny, const unsigned int nz)
bool GetPrimitive(const unsigned int i, double &a, double &b, double &c, std::vector< double > &xv, std::vector< double > &yv, std::vector< double > &zv, int &interface, double &v, double &q, double &lambda) const
void SetMinMaxNumberOfElements(const unsigned int nmin, const unsigned int nmax)
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
bool GetPlaneX(const unsigned int i, double &x, double &v) const
Retrieve the parameters of a plane at constant x.
double WeightingPotential(const double x, const double y, const double z, const std::string &label) override
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
Abstract base class for components.
Definition: Component.hh:13
std::array< bool, 3 > m_rotationSymmetric
Rotation symmetry around x-axis, y-axis, z-axis.
Definition: Component.hh:350
std::array< bool, 3 > m_mirrorPeriodic
Mirror periodicity in x, y, z.
Definition: Component.hh:346
bool m_debug
Switch on/off debugging messages.
Definition: Component.hh:341
std::array< bool, 3 > m_periodic
Simple periodicity in x, y, z.
Definition: Component.hh:344
std::string m_className
Class name.
Definition: Component.hh:329
Geometry * m_geometry
Pointer to the geometry.
Definition: Component.hh:332
bool m_ready
Ready for use?
Definition: Component.hh:338
std::array< bool, 3 > m_axiallyPeriodic
Axial periodicity in x, y, z.
Definition: Component.hh:348
virtual Solid * GetSolid(const size_t) const
Get a solid from the list.
Definition: Geometry.hh:29
virtual size_t GetNumberOfSolids() const
Return the number of solids in the geometry.
Definition: Geometry.hh:27
virtual Medium * GetMedium(const double x, const double y, const double z, const bool tesselated=false) const =0
Retrieve the medium at a given point.
Abstract base class for media.
Definition: Medium.hh:13
double GetDielectricConstant() const
Get the relative static dielectric constant.
Definition: Medium.hh:42
int GetId() const
Return the id number of the class instance.
Definition: Medium.hh:21
virtual bool IsConductor() const
Is this medium a conductor?
Definition: Medium.hh:29
bool IsDriftable() const
Is charge carrier transport enabled in this medium?
Definition: Medium.hh:74
std::string GetLabel() const
Return the label.
Definition: Solid.hh:68
bool GetCentre(double &x, double &y, double &z) const
Retrieve the centre point of the solid.
Definition: Solid.hh:71
unsigned int GetId() const
Get the ID of the solid.
Definition: Solid.hh:137
@ DielectricCharge
Definition: Solid.hh:157
@ PerpendicularField
Definition: Solid.hh:159
virtual bool SolidPanels(std::vector< Panel > &panels)=0
Retrieve the surface panels of the solid.
void Inside(const std::vector< double > &xpl, const std::vector< double > &ypl, const double x, const double y, bool &inside, bool &edge)
Definition: Polygon.cc:187
ComponentNeBem3d * gComponentNeBem3d
DoubleAc acos(const DoubleAc &f)
Definition: DoubleAc.cpp:490
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314
neBEMGLOBAL double * Lambda
Definition: neBEM.h:75
neBEMGLOBAL int * InterfaceType
Definition: neBEM.h:68
Definition: neBEM.h:155
Surface panel.
Definition: Solid.hh:11
std::vector< double > zv
Z-coordinates of vertices.
Definition: Solid.hh:19
int volume
Reference to solid to which the panel belongs.
Definition: Solid.hh:23
double a
Perpendicular vector.
Definition: Solid.hh:13
double c
Definition: Solid.hh:13
double b
Definition: Solid.hh:13
std::vector< double > xv
X-coordinates of vertices.
Definition: Solid.hh:15
std::vector< double > yv
Y-coordinates of vertices.
Definition: Solid.hh:17