Garfield++ 3.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 <map>
7#include <numeric>
8#include <vector>
9
13#include "Garfield/Polygon.hh"
14
15namespace {
16
17unsigned int NextPoint(const unsigned int i, const unsigned int n) {
18 const unsigned int j = i + 1;
19 return j < n ? j : 0;
20}
21
22unsigned int PrevPoint(const unsigned int i, const unsigned int n) {
23 return i > 0 ? i - 1 : n - 1;
24}
25
26double Mag(const std::array<double, 3>& a) {
27 return sqrt(a[0] * a[0] + a[1] * a[1] + a[2] * a[2]);
28}
29
30std::array<double, 3> UnitVector(const std::array<double, 3>& a) {
31
32 const double mag = Mag(a);
33 if (mag < 1.e-12) return a;
34 const std::array<double, 3> b = {a[0] / mag, a[1] / mag, a[2] / mag};
35 return b;
36}
37std::array<double, 3> CrossProduct(const std::array<double, 3>& u,
38 const std::array<double, 3>& v) {
39
40 const std::array<double, 3> w = {u[1] * v[2] - u[2] * v[1],
41 u[2] * v[0] - u[0] * v[2],
42 u[0] * v[1] - u[1] * v[0]};
43 return w;
44}
45
46std::array<double, 3> LocalToGlobal(const double x, const double y,
47 const double z,
48 const std::array<std::array<double, 3>, 3>& dcos,
49 const std::array<double, 3>& t) {
50
51 // Initial vector.
52 std::array<double, 4> u = {x, y, z, 1.};
53
54 std::array<std::array<double, 4>, 4> rot;
55 rot[0] = {dcos[0][0], dcos[1][0], dcos[2][0], 0.};
56 rot[1] = {dcos[0][1], dcos[1][1], dcos[2][1], 0.};
57 rot[2] = {dcos[0][2], dcos[1][2], dcos[2][2], 0.};
58 rot[3] = {0., 0., 0., 1.};
59
60 std::array<double, 4> v = {0., 0., 0., 0.};
61 for (int i = 0; i < 4; ++i) {
62 for (int j = 0; j < 4; ++j) {
63 v[i] += rot[i][j] * u[j];
64 }
65 }
66 const std::array<double, 3> a = {v[0] + t[0], v[1] + t[1], v[2] + t[2]};
67 return a;
68}
69
70/// Compute lambda for a point on a line (0 = start, 1 = end).
71double Lambda(const double x1, const double x0, const double x2,
72 const double y1, const double y0, const double y2) {
73 // Segment of zero length.
74 if ((x1 - x2) == 0. && (y1 - y2) == 0.) {
75 std::cerr << "ComponentNeBem3d::Lambda: Zero length segment.\n";
76 return 2.;
77 }
78
79 double xl = 0.;
80 const double dx1 = x0 - x1;
81 const double dy1 = y0 - y1;
82 const double dx2 = x0 - x2;
83 const double dy2 = y0 - y2;
84 if (dx1 * dx1 + dy1 * dy1 < dx2 * dx2 + dy2 * dy2) {
85 // Point nearer to (x1, y1).
86 if (fabs(y1 - y2) > fabs(x1 - x2)) {
87 xl = dy1 / (y2 - y1);
88 } else {
89 xl = dx1 / (x2 - x1);
90 }
91 } else {
92 // Point nearer to (x2, y2).
93 if (fabs(y1 - y2) > fabs(x1 - x2)) {
94 xl = 1. - dy2 / (y1 - y2);
95 } else {
96 xl = 1. - dx2 / (x1 - x2);
97 }
98 }
99 return xl;
100}
101
102/// Determine whether a point (u, v) lies on a straight line
103/// (x1, y1) to (x2, y2).
104bool OnLine(const double x1, const double y1, const double x2, const double y2,
105 const double u, const double v) {
106 // Set tolerances.
107 double epsx = 1.e-10 * std::max({fabs(x1), fabs(x2), fabs(u)});
108 double epsy = 1.e-10 * std::max({fabs(y1), fabs(y2), fabs(v)});
109 epsx = std::max(1.e-10, epsx);
110 epsy = std::max(1.e-10, epsy);
111
112 if ((fabs(x1 - u) <= epsx && fabs(y1 - v) <= epsy) ||
113 (fabs(x2 - u) <= epsx && fabs(y2 - v) <= epsy)) {
114 // Point to be examined coincides with start or end.
115 return true;
116 } else if (fabs(x1 - x2) <= epsx && fabs(y1 - y2) <= epsy) {
117 // The line (x1, y1) to (x2, y2) is in fact a point.
118 return false;
119 }
120 double xc = 0., yc = 0.;
121 if (fabs(u - x1) + fabs(v - y1) < fabs(u - x2) + fabs(v - y2)) {
122 // (u, v) is nearer to (x1, y1).
123 const double dx = (x2 - x1);
124 const double dy = (y2 - y1);
125 const double xl = ((u - x1) * dx + (v - y1) * dy) / (dx * dx + dy * dy);
126 if (xl < 0.) {
127 xc = x1;
128 yc = y1;
129 } else if (xl > 1.) {
130 xc = x2;
131 yc = y2;
132 } else {
133 xc = x1 + xl * dx;
134 yc = y1 + xl * dy;
135 }
136 } else {
137 // (u, v) is nearer to (x2, y2).
138 const double dx = (x1 - x2);
139 const double dy = (y1 - y2);
140 const double xl = ((u - x2) * dx + (v - y2) * dy) / (dx * dx + dy * dy);
141 if (xl < 0.) {
142 xc = x2;
143 yc = y2;
144 } else if (xl > 1.) {
145 xc = x1;
146 yc = y1;
147 } else {
148 xc = x2 + xl * dx;
149 yc = y2 + xl * dy;
150 }
151 }
152 // See whether the point is on the line.
153 if (fabs(u - xc) < epsx && fabs(v - yc) < epsy) {
154 return true;
155 }
156 return false;
157}
158
159/// Determine whether the 2 straight lines (x1, y1) to (x2, y2)
160/// and (u1, v1) to (u2, v2) cross at an intermediate point for both lines.
161bool Crossing(const double x1, const double y1, const double x2,
162 const double y2, const double u1, const double v1,
163 const double u2, const double v2, double& xc, double& yc) {
164 /// Matrix to compute the crossing point.
165 std::array<std::array<double, 2>, 2> a;
166 a[0][0] = y2 - y1;
167 a[0][1] = v2 - v1;
168 a[1][0] = x1 - x2;
169 a[1][1] = u1 - u2;
170 const double det = a[0][0] * a[1][1] - a[1][0] * a[0][1];
171 // Initial values.
172 xc = 0.;
173 yc = 0.;
174 // Set tolerances.
175 double epsx = 1.e-10 * std::max({fabs(x1), fabs(x2), fabs(u1), fabs(u2)});
176 double epsy = 1.e-10 * std::max({fabs(y1), fabs(y2), fabs(v1), fabs(v2)});
177 epsx = std::max(epsx, 1.e-10);
178 epsy = std::max(epsy, 1.e-10);
179 // Check for a point of one line located on the other line.
180 if (OnLine(x1, y1, x2, y2, u1, v1)) {
181 xc = u1;
182 yc = v1;
183 return true;
184 } else if (OnLine(x1, y1, x2, y2, u2, v2)) {
185 xc = u2;
186 yc = v2;
187 return true;
188 } else if (OnLine(u1, v1, u2, v2, x1, y1)) {
189 xc = x1;
190 yc = y1;
191 return true;
192 } else if (OnLine(u1, v1, u2, v2, x2, y2)) {
193 xc = x2;
194 yc = y2;
195 return true;
196 } else if (fabs(det) < epsx * epsy) {
197 // Parallel, non-touching.
198 return false;
199 }
200 // Crossing, non-trivial lines: solve crossing equations.
201 const double aux = a[1][1];
202 a[1][1] = a[0][0] / det;
203 a[0][0] = aux / det;
204 a[1][0] = -a[1][0] / det;
205 a[0][1] = -a[0][1] / det;
206 // Compute crossing point.
207 xc = a[0][0] * (x1 * y2 - x2 * y1) + a[1][0] * (u1 * v2 - u2 * v1);
208 yc = a[0][1] * (x1 * y2 - x2 * y1) + a[1][1] * (u1 * v2 - u2 * v1);
209 // See whether the crossing point is on both lines.
210 if (OnLine(x1, y1, x2, y2, xc, yc) && OnLine(u1, v1, u2, v2, xc, yc)) {
211 // Intersecting lines.
212 return true;
213 }
214 // Crossing point not on both lines.
215 return false;
216}
217
218void AddPoints(const std::vector<double>& xp1, const std::vector<double>& yp1,
219 const std::vector<double>& xp2, const std::vector<double>& yp2,
220 std::vector<double>& xl, std::vector<double>& yl,
221 std::vector<int>& flags, std::vector<double>& qs,
222 const double epsx, const double epsy) {
223 struct Point {
224 double x;
225 double y;
226 int flag;
227 double q;
228 };
229
230 std::vector<Point> points;
231
232 const unsigned int np1 = xp1.size();
233 const unsigned int np2 = xp2.size();
234 for (unsigned int i = 0; i < np1; ++i) {
235 const double xi0 = xp1[i];
236 const double yi0 = yp1[i];
237 const double xi1 = xp1[NextPoint(i, np1)];
238 const double yi1 = yp1[NextPoint(i, np1)];
239 // Add the vertex.
240 Point p1;
241 p1.x = xi0;
242 p1.y = yi0;
243 p1.flag = 1;
244 p1.q = 0.;
245 // If also on 2 or vertex of 2, flag it as crossing or foreign.
246 for (unsigned int j = 0; j < np2; ++j) {
247 const double xj0 = xp2[j];
248 const double yj0 = yp2[j];
249 if (fabs(xj0 - xi0) < epsx && fabs(yj0 - yi0) < epsy) {
250 p1.flag = 2;
251 }
252 const double xj1 = xp2[NextPoint(j, np2)];
253 const double yj1 = yp2[NextPoint(j, np2)];
254 if (OnLine(xj0, yj0, xj1, yj1, xi0, yi0) &&
255 (fabs(xj0 - xi0) > epsx || fabs(yj0 - yi0) > epsy) &&
256 (fabs(xj1 - xi0) > epsx || fabs(yj1 - yi0) > epsy)) {
257 p1.flag = 3;
258 }
259 }
260 points.push_back(std::move(p1));
261 // Go over the line segments of the other polygon.
262 std::vector<Point> pointsOther;
263 for (unsigned int j = 0; j < np2; ++j) {
264 const double xj0 = xp2[j];
265 const double yj0 = yp2[j];
266 // Add vertices of 2 that are on this line.
267 if (OnLine(xi0, yi0, xi1, yi1, xj0, yj0) &&
268 (fabs(xi0 - xj0) > epsx || fabs(yi0 - yj0) > epsy) &&
269 (fabs(xi1 - xj0) > epsx || fabs(yi1 - yj0) > epsy)) {
270 Point p2;
271 p2.x = xj0;
272 p2.y = yj0;
273 p2.flag = 2;
274 pointsOther.push_back(std::move(p2));
275 }
276 const double xj1 = xp2[NextPoint(j, np2)];
277 const double yj1 = yp2[NextPoint(j, np2)];
278 // Add crossing points.
279 double xc = 0., yc = 0.;
280 bool add = Crossing(xi0, yi0, xi1, yi1, xj0, yj0, xj1, yj1, xc, yc);
281 if (add) {
282 if ((fabs(xi0 - xc) < epsx && fabs(yi0 - yc) < epsy) ||
283 (fabs(xi1 - xc) < epsx && fabs(yi1 - yc) < epsy) ||
284 (fabs(xj0 - xc) < epsx && fabs(yj0 - yc) < epsy) ||
285 (fabs(xj1 - xc) < epsx && fabs(yj1 - yc) < epsy)) {
286 add = false;
287 }
288 if ((fabs(xi0 - xj0) < epsx && fabs(yi0 - yj0) < epsy) ||
289 (fabs(xi0 - xj1) < epsx && fabs(yi0 - yj1) < epsy) ||
290 (fabs(xi1 - xj0) < epsx && fabs(yi1 - yj0) < epsy) ||
291 (fabs(xi1 - xj1) < epsx && fabs(yi1 - yj1) < epsy)) {
292 add = false;
293 }
294 }
295 if (add) {
296 Point p2;
297 p2.x = xc;
298 p2.y = yc;
299 p2.flag = 3;
300 pointsOther.push_back(std::move(p2));
301 }
302 }
303 // Compute the lambdas for these points.
304 for (auto& p : pointsOther) {
305 p.q = Lambda(xi0, p.x, xi1, yi0, p.y, yi1);
306 }
307 // Sort the list by using the lambdas.
308 std::sort(
309 pointsOther.begin(), pointsOther.end(),
310 [](const Point& lhs, const Point& rhs) { return (lhs.q < rhs.q); });
311 points.insert(points.end(), pointsOther.begin(), pointsOther.end());
312 }
313
314 for (const Point& p : points) {
315 xl.push_back(p.x);
316 yl.push_back(p.y);
317 flags.push_back(p.flag);
318 qs.push_back(p.q);
319 }
320}
321
322/// Determine whether 2 panels are equal.
323bool Equal(const Garfield::Panel& panel1, const Garfield::Panel& panel2,
324 const double epsx, const double epsy) {
325 const auto& xp1 = panel1.xv;
326 const auto& yp1 = panel1.yv;
327 const auto& xp2 = panel2.xv;
328 const auto& yp2 = panel2.yv;
329 if (xp1.empty() || xp2.empty()) return false;
330 const unsigned int np1 = xp1.size();
331 const unsigned int np2 = xp2.size();
332
333 // Compare all points of 1 with all points of 2.
334 for (unsigned int i = 0; i < np1; ++i) {
335 // Loop over 2 until a match is found.
336 bool match = false;
337 for (unsigned int j = 0; j < np2; ++j) {
338 if (fabs(xp2[j] - xp1[i]) < epsx && fabs(yp2[j] - yp1[i]) < epsy) {
339 match = true;
340 break;
341 }
342 const unsigned int jj = NextPoint(j, np2);
343 if (OnLine(xp2[j], yp2[j], xp2[jj], yp2[jj], xp1[i], yp1[i])) {
344 match = true;
345 break;
346 }
347 }
348 if (!match) return false;
349 }
350
351 // Compare all points of 2 with all points of 1.
352 for (unsigned int i = 0; i < np2; ++i) {
353 // Loop over 1 until a match is found.
354 bool match = false;
355 for (unsigned int j = 0; j < np1; ++j) {
356 if (fabs(xp2[i] - xp1[j]) < epsx && fabs(yp2[i] - yp1[j]) < epsy) {
357 match = true;
358 break;
359 }
360 const unsigned int jj = NextPoint(j, np1);
361 if (OnLine(xp1[j], yp1[j], xp1[jj], yp1[jj], xp2[i], yp2[i])) {
362 match = true;
363 break;
364 }
365 }
366 if (!match) return false;
367 }
368
369 // If we get this far, the curves are the same.
370 return true;
371}
372
373}
374
375namespace Garfield {
376
378 m_className = "ComponentNeBem3d";
379}
380
381void ComponentNeBem3d::ElectricField(const double x, const double y,
382 const double z, double& ex, double& ey,
383 double& ez, double& v, Medium*& m,
384 int& status) {
385 ex = ey = ez = v = 0.;
386 status = 0;
387 // Check if the requested point is inside a medium
388 m = GetMedium(x, y, z);
389 if (!m) {
390 status = -6;
391 return;
392 }
393
394 if (!m_ready) {
395 if (!Initialise()) {
396 std::cerr << m_className << "::ElectricField: Initialisation failed.\n";
397 status = -11;
398 return;
399 }
400 m_ready = true;
401 }
402}
403
404void ComponentNeBem3d::ElectricField(const double x, const double y,
405 const double z, double& ex, double& ey,
406 double& ez, Medium*& m, int& status) {
407 double v = 0.;
408 ElectricField(x, y, z, ex, ey, ez, v, m, status);
409}
410
411bool ComponentNeBem3d::GetVoltageRange(double& vmin, double& vmax) {
412 // Voltage and other bc have to come from the solids
413 vmin = vmax = 0;
414 return true;
415}
416
417void ComponentNeBem3d::SetTargetElementSize(const double length) {
418
419 if (length < MinDist) {
420 std::cerr << m_className << "::SetTargetElementSize: Value must be > "
421 << MinDist << ".\n";
422 return;
423 }
424 m_targetElementSize = length;
425}
426
428 // Reset the lists.
429 m_primitives.clear();
430 m_elements.clear();
431
432 if (!m_geometry) {
433 std::cerr << m_className << "::Initialise: Geometry not set.\n";
434 return false;
435 }
436 // Be sure we won't have intersections with the bounding box.
437 // TODO! Loop over the solids and call PLACYE, PLACHE, PLABXE, PLASPE, ...
438 // Loop over the solids.
439 std::map<int, Solid*> solids;
440 std::map<int, Solid::BoundaryCondition> bc;
441 std::map<int, double> volt;
442 std::map<int, double> eps;
443 std::map<int, double> charge;
444 const unsigned int nSolids = m_geometry->GetNumberOfSolids();
445 std::vector<Panel> panelsIn;
446 for (unsigned int i = 0; i < nSolids; ++i) {
447 Medium* medium = nullptr;
448 const auto solid = m_geometry->GetSolid(i, medium);
449 if (!solid) continue;
450 // Get the panels.
451 solid->SolidPanels(panelsIn);
452 // Get the boundary condition.
453 const auto id = solid->GetId();
454 solids[id] = solid;
455 bc[id] = solid->GetBoundaryConditionType();
456 volt[id] = solid->GetBoundaryPotential();
457 charge[id] = solid->GetBoundaryChargeDensity();
458 if (!medium) {
459 eps[id] = 1.;
460 } else {
461 eps[id] = medium->GetDielectricConstant();
462 }
463 }
464 // Apply cuts.
465 // CALL CELSCT('APPLY')
466 // Reduce to basic periodic copy.
467 // CALL BEMBAS
468
469 // Find contact panels and split into primitives.
470
471 // *---------------------------------------------------------------------
472 // * PLABEM - Prepares panels for BEM applications: removes the contacts
473 // * and cuts polygons to rectangles and right-angle triangles.
474 // *---------------------------------------------------------------------
475
476 // Establish tolerances.
477 const double epsang = 1.e-6; // BEMEPA
478 const double epsxyz = 1.e-6; // BEMEPD
479 // CALL EPSSET('SET',EPSXYZ,EPSXYZ,EPSXYZ)
480
481 const unsigned int nIn = panelsIn.size();
482 if (m_debug) {
483 std::cout << m_className << "::Initialise: Retrieved "
484 << nIn << " panels from the solids.\n";
485 }
486 // Keep track of which panels have been processed.
487 std::vector<bool> mark(nIn, false);
488 // Pick up panels which coincide potentially.
489 for (unsigned int i = 0; i < nIn; ++i) {
490 // Skip panels already done.
491 if (mark[i]) continue;
492 // Fetch panel parameters.
493 const double a1 = panelsIn[i].a;
494 const double b1 = panelsIn[i].b;
495 const double c1 = panelsIn[i].c;
496 const auto& xp1 = panelsIn[i].xv;
497 const auto& yp1 = panelsIn[i].yv;
498 const auto& zp1 = panelsIn[i].zv;
499 const unsigned int np1 = xp1.size();
500 // Establish its norm and offset.
501 const double d1 = a1 * xp1[0] + b1 * yp1[0] + c1 * zp1[0];
502 if (m_debug) {
503 std::cout << " Panel " << i << "\n Norm vector: "
504 << a1 << ", " << b1 << ", " << c1 << ", " << d1 << ".\n";
505 }
506 // Rotation matrix.
507 std::array<std::array<double, 3>, 3> rot;
508 if (fabs(c1) <= fabs(a1) && fabs(c1) <= fabs(b1)) {
509 // Rotation: removing C
510 rot[0][0] = b1 / sqrt(a1 * a1 + b1 * b1);
511 rot[0][1] = -a1 / sqrt(a1 * a1 + b1 * b1);
512 rot[0][2] = 0.;
513 } else if (fabs(b1) <= fabs(a1) && fabs(b1) <= fabs(c1)) {
514 // Rotation: removing B
515 rot[0][0] = c1 / sqrt(a1 * a1 + c1 * c1);
516 rot[0][1] = 0.;
517 rot[0][2] = -a1 / sqrt(a1 * a1 + c1 * c1);
518 } else {
519 // Rotation: removing A
520 rot[0][0] = 0.;
521 rot[0][1] = c1 / sqrt(b1 * b1 + c1 * c1);
522 rot[0][2] = -b1 / sqrt(b1 * b1 + c1 * c1);
523 }
524 rot[2][0] = a1;
525 rot[2][1] = b1;
526 rot[2][2] = c1;
527 rot[1][0] = rot[2][1] * rot[0][2] - rot[2][2] * rot[0][1];
528 rot[1][1] = rot[2][2] * rot[0][0] - rot[2][0] * rot[0][2];
529 rot[1][2] = rot[2][0] * rot[0][1] - rot[2][1] * rot[0][0];
530 // Rotate to the x, y plane.
531 std::vector<double> xp(np1, 0.);
532 std::vector<double> yp(np1, 0.);
533 std::vector<double> zp(np1, 0.);
534 double zm = 0.;
535 for (unsigned int k = 0; k < np1; ++k) {
536 xp[k] = rot[0][0] * xp1[k] + rot[0][1] * yp1[k] + rot[0][2] * zp1[k];
537 yp[k] = rot[1][0] * xp1[k] + rot[1][1] * yp1[k] + rot[1][2] * zp1[k];
538 zp[k] = rot[2][0] * xp1[k] + rot[2][1] * yp1[k] + rot[2][2] * zp1[k];
539 zm += zp[k];
540 }
541 zm /= np1;
542 // Store it.
543 std::vector<Panel> newPanels;
544 std::vector<int> vol1;
545 std::vector<int> vol2;
546 Panel panel1 = panelsIn[i];
547 panel1.xv = xp;
548 panel1.yv = yp;
549 panel1.zv = zp;
550 vol1.push_back(panel1.volume);
551 vol2.push_back(-1);
552 newPanels.push_back(std::move(panel1));
553 // Pick up all matching planes.
554 for (unsigned int j = i + 1; j < nIn; ++j) {
555 if (mark[j]) continue;
556 const double a2 = panelsIn[j].a;
557 const double b2 = panelsIn[j].b;
558 const double c2 = panelsIn[j].c;
559 const auto& xp2 = panelsIn[j].xv;
560 const auto& yp2 = panelsIn[j].yv;
561 const auto& zp2 = panelsIn[j].zv;
562 const unsigned int np2 = xp2.size();
563 // See whether this matches the first.
564 const double d2 = a2 * xp2[0] + b2 * yp2[0] + c2 * zp2[0];
565 // Inner product.
566 const double dot = a1 * a2 + b1 * b2 + c1 * c2;
567 // Offset between the two planes.
568 const double offset = d1 - d2 * dot;
569 if (fabs(fabs(dot) - 1.) > epsang || fabs(offset) > epsxyz) continue;
570 // Found a match.
571 mark[j] = true;
572 if (m_debug) std::cout << " Match with panel " << j << ".\n";
573 // Rotate this plane too.
574 xp.assign(np2, 0.);
575 yp.assign(np2, 0.);
576 zp.assign(np2, 0.);
577 zm = 0.;
578 for (unsigned int k = 0; k < np2; ++k) {
579 xp[k] = rot[0][0] * xp2[k] + rot[0][1] * yp2[k] + rot[0][2] * zp2[k];
580 yp[k] = rot[1][0] * xp2[k] + rot[1][1] * yp2[k] + rot[1][2] * zp2[k];
581 zp[k] = rot[2][0] * xp2[k] + rot[2][1] * yp2[k] + rot[2][2] * zp2[k];
582 zm += zp[k];
583 }
584 zm /= np2;
585 // Store it.
586 Panel panel2 = panelsIn[j];
587 panel2.xv = xp;
588 panel2.yv = yp;
589 panel2.zv = zp;
590 vol1.push_back(panel2.volume);
591 vol2.push_back(-1);
592 newPanels.push_back(std::move(panel2));
593 }
594 std::vector<bool> obsolete(newPanels.size(), false);
595 // Cut them as long as needed till no contacts remain.
596 unsigned int jmin = 0;
597 bool change = true;
598 while (change) {
599 change = false;
600 const unsigned int n = newPanels.size();
601 for (unsigned int j = 0; j < n; ++j) {
602 if (obsolete[j] || j < jmin) continue;
603 if (vol1[j] >= 0 && vol2[j] >= 0) continue;
604 const auto& panelj = newPanels[j];
605 for (unsigned int k = j + 1; k < n; ++k) {
606 if (obsolete[k]) continue;
607 if (vol1[k] >= 0 && vol2[k] >= 0) continue;
608 const auto& panelk = newPanels[k];
609 if (m_debug) std::cout << " Cutting " << j << ", " << k << ".\n";
610 // Separate contact and non-contact areas.
611 std::vector<Panel> panelsOut;
612 std::vector<int> itypo;
613 EliminateOverlaps(panelj, panelk, panelsOut, itypo);
614 const unsigned int nOut = panelsOut.size();
615 if (nOut == 2) {
616 // TODO: retrieve epsx, epsy from overlap finding?
617 const double epsx = epsxyz;
618 const double epsy = epsxyz;
619 // If there are just 2 panels, see whether there is a new one.
620 const bool equal1 = Equal(panelj, panelsOut[0], epsx, epsy);
621 const bool equal2 = Equal(panelj, panelsOut[1], epsx, epsy);
622 const bool equal3 = Equal(panelk, panelsOut[0], epsx, epsy);
623 const bool equal4 = Equal(panelk, panelsOut[1], epsx, epsy);
624 if ((equal1 || equal3) && (equal2 || equal4)) {
625 if (m_debug) {
626 std::cout << " Original and new panels are identical.\n";
627 }
628 } else {
629 change = true;
630 }
631 } else {
632 change = true;
633 }
634 if (m_debug) std::cout << " Change flag: " << change << ".\n";
635 // If there is no change, keep the two panels and proceed.
636 if (!change) continue;
637 // Flag the existing panels as inactive.
638 obsolete[j] = true;
639 obsolete[k] = true;
640
641 // Add the new panels.
642 for (unsigned int l = 0; l < nOut; ++l) {
643 if (itypo[l] == 1) {
644 vol1.push_back(std::max(vol1[j], vol2[j]));
645 vol2.push_back(-1);
646 } else if (itypo[l] == 2) {
647 vol1.push_back(std::max(vol1[k], vol2[k]));
648 vol2.push_back(-1);
649 } else {
650 vol1.push_back(std::max(vol1[j], vol2[j]));
651 vol2.push_back(std::max(vol1[k], vol2[k]));
652 }
653 newPanels.push_back(std::move(panelsOut[l]));
654 obsolete.push_back(false);
655 }
656 jmin = j + 1;
657 // Restart the loops.
658 break;
659 }
660 if (change) break;
661 }
662 }
663 // And rotate the panels back in place.
664 const unsigned int nNew = newPanels.size();
665 for (unsigned int j = 0; j < nNew; ++j) {
666 if (obsolete[j]) continue;
667 // Examine the boundary conditions.
668 int interfaceType = 0;
669 double potential = 0.;
670 double lambda = 0.;
671 double chargeDensity = 0.;
672 if (m_debug) {
673 std::cout << " Volume 1: " << vol1[j]
674 << ". Volume 2: " << vol2[j] << ".\n";
675 }
676 if (vol1[j] < 0 && vol2[j] < 0) {
677 // Shouldn't happen.
678 continue;
679 } else if (vol1[j] < 0 || vol2[j] < 0) {
680 // Interface between a solid and vacuum/background.
681 const auto vol = vol1[j] < 0 ? vol2[j] : vol1[j];
682 interfaceType = InterfaceType(bc[vol]);
683 if (bc[vol] == Solid::Dielectric ||
684 bc[vol] == Solid::DielectricCharge) {
685 if (fabs(eps[vol] - 1.) < 1.e-6) {
686 // Same epsilon on both sides. Skip.
687 interfaceType = 0;
688 } else {
689 lambda = (eps[vol] - 1.) / (eps[vol] + 1.);
690 }
691 } else if (bc[vol] == Solid::Voltage) {
692 potential = volt[vol];
693 }
694 if (bc[vol] == Solid::Charge || bc[vol] == Solid::DielectricCharge) {
695 chargeDensity = charge[vol];
696 }
697 } else {
698 const auto bc1 = bc[vol1[j]];
699 const auto bc2 = bc[vol2[j]];
700 if (bc1 == Solid::Voltage || bc1 == Solid::Charge ||
701 bc1 == Solid::Float) {
702 interfaceType = InterfaceType(bc1);
703 // First volume is a conductor. Other volume must be a dielectric.
704 if (bc2 == Solid::Dielectric || bc2 == Solid::DielectricCharge) {
705 if (bc1 == Solid::Voltage) {
706 potential = volt[vol1[j]];
707 } else if (bc1 == Solid::Charge) {
708 chargeDensity = charge[vol1[j]];
709 }
710 } else {
711 interfaceType = -1;
712 }
713 if (bc1 == Solid::Voltage && bc2 == Solid::Voltage) {
714 const double v1 = volt[vol1[j]];
715 const double v2 = volt[vol2[j]];
716 if (fabs(v1 - v2) < 1.e-6 * (1. + fabs(v1) + fabs(v2))) {
717 interfaceType = 0;
718 }
719 }
720 } else if (bc1 == Solid::Dielectric ||
722 interfaceType = InterfaceType(bc2);
723 // First volume is a dielectric.
724 if (bc2 == Solid::Voltage) {
725 potential = volt[vol2[j]];
726 } else if (bc2 == Solid::Charge) {
727 chargeDensity = charge[vol2[j]];
728 } else if (bc2 == Solid::Dielectric ||
730 const double eps1 = eps[vol1[j]];
731 const double eps2 = eps[vol2[j]];
732 if (fabs(eps1 - eps2) < 1.e-6 * (1. + fabs(eps1) + fabs(eps2))) {
733 // Same epsilon. Skip.
734 interfaceType = 0;
735 } else {
736 lambda = (eps1 - eps2) / (eps1 + eps2);
737 }
738 }
739 }
740 }
741 if (m_debug) {
742 if (interfaceType < 0) {
743 std::cout << " Conflicting boundary conditions. Skip.\n";
744 } else if (interfaceType < 1) {
745 std::cout << " Trivial interface. Skip.\n";
746 } else if (interfaceType > 5) {
747 std::cout << " Interface type " << interfaceType
748 << " is not implemented. Skip.\n";
749 }
750 }
751 if (interfaceType < 1 || interfaceType > 5) continue;
752
753 std::vector<Panel> panelsOut;
754 // Reduce to rectangles and right-angle triangles.
755 if (m_debug) std::cout << " Creating primitives.\n";
756 MakePrimitives(newPanels[j], panelsOut);
757 // Loop over the rectangles and triangles.
758 for (auto& panel : panelsOut) {
759 const auto& up = panel.xv;
760 const auto& vp = panel.yv;
761 const auto& wp = panel.zv;
762 const unsigned int np = up.size();
763 // Rotate.
764 xp.assign(np, 0.);
765 yp.assign(np, 0.);
766 zp.assign(np, 0.);
767 for (unsigned int k = 0; k < np; ++k) {
768 xp[k] = rot[0][0] * up[k] + rot[1][0] * vp[k] + rot[2][0] * wp[k];
769 yp[k] = rot[0][1] * up[k] + rot[1][1] * vp[k] + rot[2][1] * wp[k];
770 zp[k] = rot[0][2] * up[k] + rot[1][2] * vp[k] + rot[2][2] * wp[k];
771 }
772 Primitive primitive;
773 primitive.a = panel.a;
774 primitive.b = panel.b;
775 primitive.c = panel.c;
776 primitive.xv = xp;
777 primitive.yv = yp;
778 primitive.zv = zp;
779 primitive.v = potential;
780 primitive.q = chargeDensity;
781 primitive.lambda = lambda;
782 primitive.interface = interfaceType;
783 // Set the requested discretization level (target element size).
784 primitive.elementSize = -1.;
785 if (solids.find(vol1[j]) != solids.end()) {
786 const auto solid = solids[vol1[j]];
787 if (solid) {
788 primitive.elementSize = solid->GetDiscretisationLevel(panel);
789 }
790 }
791 m_primitives.push_back(std::move(primitive));
792 }
793 }
794 }
795
796 if (m_debug) {
797 std::cout << m_className << "::Initialise:\n"
798 << " Created " << m_primitives.size() << " primitives.\n";
799 }
800
801 // Discretize the primitives.
802 for (const auto& primitive : m_primitives) {
803 const auto nVertices = primitive.xv.size();
804 if (nVertices < 2 || nVertices > 4) continue;
805 std::vector<Element> elements;
806 // Get the target element size.
807 double targetSize = primitive.elementSize;
808 if (targetSize < MinDist) targetSize = m_targetElementSize;
809 if (nVertices == 2) {
810 DiscretizeWire(primitive, targetSize, elements);
811 } else if (nVertices == 3) {
812 DiscretizeTriangle(primitive, targetSize, elements);
813 } else if (nVertices == 4) {
814 DiscretizeRectangle(primitive, targetSize, elements);
815 }
816 for (auto& element: elements) {
817 element.interface = primitive.interface;
818 element.lambda = primitive.lambda;
819 element.bc = primitive.v;
820 element.assigned = primitive.q;
821 element.solution = 0.;
822 }
823 m_elements.insert(m_elements.end(),
824 std::make_move_iterator(elements.begin()),
825 std::make_move_iterator(elements.end()));
826 }
827 return true;
828}
829
830int ComponentNeBem3d::InterfaceType(const Solid::BoundaryCondition bc) const {
831
832 // 0: To be skipped
833 // 1: Conductor-dielectric
834 // 2: Conductor with known charge
835 // 3: Conductor at floating potential
836 // 4: Dielectric-dielectric
837 // 5: Dielectric with given surface charge
838 //
839 // Check dielectric-dielectric formulation in
840 // Jaydeep P. Bardhan,
841 // Numerical solution of boundary-integral equations for molecular
842 // electrostatics, J. Chem. Phys. 130, 094102 (2009)
843 // https://doi.org/10.1063/1.3080769
844
845 switch (bc) {
846 case Solid::Voltage:
847 return 1;
848 break;
849 case Solid::Charge:
850 return 2;
851 break;
852 case Solid::Float:
853 return 3;
854 break;
856 return 4;
857 break;
859 return 5;
860 break;
862 // Symmetry boundary, E parallel.
863 return 6;
864 break;
866 // Symmetry boundary, E perpendicular.
867 return 7;
868 break;
869 default:
870 break;
871 }
872 return 0;
873}
874
875bool ComponentNeBem3d::DiscretizeWire(const Primitive& primitive,
876 const double targetSize, std::vector<Element>& elements) const {
877
878 const double dx = primitive.xv[1] - primitive.xv[0];
879 const double dy = primitive.yv[1] - primitive.yv[0];
880 const double dz = primitive.zv[1] - primitive.zv[0];
881 const double lw = sqrt(dx * dx + dy * dy + dz * dz);
882 // Determine the number of segments along the wire.
883 unsigned int nSegments = NbOfSegments(lw, targetSize);
884 const double elementSize = lw / nSegments;
885
886 // Determine the direction cosines.
887 // The z axis of the local coordinate system is along the wire.
888 const std::array<double, 3> zu = {dx / lw, dy / lw, dz / lw};
889 // Make the x axis orthogonal in the two largest components.
890 std::array<double, 3> xu;
891 if (fabs(zu[0]) >= fabs(zu[2]) && fabs(zu[1]) >= fabs(zu[2])) {
892 xu = {-zu[1], zu[0], 0.};
893 } else if (fabs(zu[0]) >= fabs(zu[1]) && fabs(zu[2]) >= fabs(zu[1])) {
894 xu = {-zu[2], 0., zu[0]};
895 } else {
896 xu = {0., zu[2], -zu[1]};
897 }
898 xu = UnitVector(xu);
899 // The y axis is given by the vectorial product of z and x.
900 const std::array<double, 3> yu = UnitVector(CrossProduct(zu, xu));
901 const std::array<std::array<double, 3>, 3> dcos = {xu, yu, zu};
902
903 const double xincr = dx / nSegments;
904 const double yincr = dy / nSegments;
905 const double zincr = dz / nSegments;
906
907 // TODO!
908 const double radius = 1.;
909 const double dA = TwoPi * radius * elementSize;
910 for (unsigned int i = 0; i < nSegments; ++i) {
911 const double x0 = primitive.xv[0] + i * xincr;
912 const double y0 = primitive.yv[0] + i * yincr;
913 const double z0 = primitive.zv[0] + i * zincr;
914 Element element;
915 element.origin = {x0 + 0.5 * xincr, y0 + 0.5 * yincr, z0 + 0.5 * zincr};
916 element.xv = {x0, x0 + xincr};
917 element.yv = {y0, y0 + yincr};
918 element.zv = {z0, z0 + zincr};
919 element.lx = radius;
920 element.lz = elementSize;
921 element.dA = dA;
922 element.dcos = dcos;
923 // Modify to be on surface?
924 element.collocationPoint = element.origin;
925 elements.push_back(std::move(element));
926 }
927 return true;
928}
929
930bool ComponentNeBem3d::DiscretizeTriangle(const Primitive& primitive,
931 const double targetSize, std::vector<Element>& elements) const {
932
933 // Origin of the local coordinate system is at the right angle corner.
934 std::array<double, 3> corner = {primitive.xv[1], primitive.yv[1],
935 primitive.zv[1]};
936 // Determine the direction cosines.
937 const double dx1 = primitive.xv[0] - primitive.xv[1];
938 const double dy1 = primitive.yv[0] - primitive.yv[1];
939 const double dz1 = primitive.zv[0] - primitive.zv[1];
940 const double dx2 = primitive.xv[2] - primitive.xv[1];
941 const double dy2 = primitive.yv[2] - primitive.yv[1];
942 const double dz2 = primitive.zv[2] - primitive.zv[1];
943 // Normal vector of the primitive.
944 std::array<double, 3> nu = {primitive.a, primitive.b, primitive.c};
945 nu = UnitVector(nu);
946 // int flagDC = 1;
947 // We begin with trial 1: one of the possible orientations.
948 double lx = sqrt(dx1 * dx1 + dy1 * dy1 + dz1 * dz1);
949 double lz = sqrt(dx2 * dx2 + dy2 * dy2 + dz2 * dz2);
950 std::array<double, 3> xu = {dx1 / lx, dy1 / lx, dz1 / lx};
951 std::array<double, 3> zu = {dx2 / lz, dy2 / lz, dz2 / lz};
952 std::array<double, 3> yu = CrossProduct(zu, xu);
953 constexpr double tol = 1.e-3;
954 if ((fabs(yu[0] - nu[0]) > tol) || (fabs(yu[1] - nu[1]) > tol) ||
955 (fabs(yu[2] - nu[2]) > tol)) {
956 // flagDC = 2;
957 // Try the other orientation.
958 std::swap(lx, lz);
959 xu = {dx2 / lx, dy2 / lx, dz2 / lx};
960 zu = {dx1 / lz, dy1 / lz, dz1 / lz};
961 yu = CrossProduct(zu, xu);
962 if ((fabs(yu[0] - nu[0]) > tol) || (fabs(yu[1] - nu[1]) > tol) ||
963 (fabs(yu[2] - nu[2]) > tol)) {
964 // No other possibility, search failed.
965 std::cerr << m_className << "::DiscretizeTriangle:\n"
966 << " Could not establish direction vectors.\n";
967 return false;
968 }
969 }
970 std::array<std::array<double, 3>, 3> dcos = {xu, yu, zu};
971
972 // TODO!
973 /*
974 double eps1 = primitive.eps1;
975 double eps2 = primitive.eps2;
976 if (flagDC == 1) std::swap(eps1, eps2);
977 */
978
979 // Determine the number of elements.
980 unsigned int nx = NbOfSegments(lx, targetSize);
981 unsigned int nz = NbOfSegments(lz, targetSize);
982 double elementSizeX = lx / nx;
983 double elementSizeZ = lz / nz;
984
985 // Analyse the aspect ratio.
986 double ar = elementSizeX / elementSizeZ;
987 if (ar > 10.) {
988 // Reduce nz.
989 elementSizeZ = 0.1 * elementSizeX;
990 nz = std::max(static_cast<int>(lz / elementSizeZ), 1);
991 elementSizeZ = lz / nz;
992 }
993 if (ar < 0.1) {
994 // Reduce nx.
995 elementSizeX = 0.1 * elementSizeZ;
996 nx = std::max(static_cast<int>(lx / elementSizeX), 1);
997 elementSizeX = lx / nx;
998 }
999 const double dxdz = lx / lz;
1000 for (unsigned int k = 0; k < nz; ++k) {
1001 // Consider the k-th row.
1002 const double zlo = k * elementSizeZ;
1003 const double zhi = (k + 1) * elementSizeZ;
1004 double xlo = (lz - zlo) * dxdz;
1005 double xhi = (lz - zhi) * dxdz;
1006 // Triangular element on the k-th row.
1007 Element triangle;
1008 triangle.origin = LocalToGlobal(xhi, 0., zlo, dcos, corner);
1009 // Assign element values.
1010 const double triangleSizeX = xlo - xhi;
1011 triangle.lx = triangleSizeX;
1012 triangle.lz = elementSizeZ;
1013 triangle.dA = 0.5 * triangleSizeX * elementSizeZ;
1014 triangle.dcos = dcos;
1015 // Calculate the element vertices.
1016 const double xv0 = triangle.origin[0];
1017 const double yv0 = triangle.origin[1];
1018 const double zv0 = triangle.origin[2];
1019 const double xv1 = xv0 + triangleSizeX * xu[0];
1020 const double yv1 = yv0 + triangleSizeX * xu[1];
1021 const double zv1 = zv0 + triangleSizeX * xu[2];
1022 const double xv2 = xv0 + elementSizeZ * zu[0];
1023 const double yv2 = yv0 + elementSizeZ * zu[1];
1024 const double zv2 = zv0 + elementSizeZ * zu[2];
1025 // Assign vertices of the element.
1026 triangle.xv = {xv0, xv1, xv2};
1027 triangle.yv = {yv0, yv1, yv2};
1028 triangle.zv = {zv0, zv1, zv2};
1029
1030 // Boundary condition is applied at the barycentre, not at the origin
1031 // of the element coordinate system which is at the right-angled corner.
1032 const double xb = triangleSizeX / 3.;
1033 const double yb = 0.;
1034 const double zb = elementSizeZ / 3.;
1035 triangle.collocationPoint = LocalToGlobal(xb, yb, zb, dcos, triangle.origin);
1036 elements.push_back(std::move(triangle));
1037 // No rectangular element on the last row.
1038 if (k == nz - 1) continue;
1039
1040 // Determine number and size in x of the rectangular elements.
1041 const int nRect = xhi <= elementSizeX ? 1 : (int)(xhi / elementSizeX);
1042 const double rectSizeX = xhi / nRect;
1043 const double zc = 0.5 * (zlo + zhi);
1044 for (int i = 0; i < nRect; ++i) {
1045 // Centroid of the rectangular element.
1046 const double xc = (i + 0.5) * rectSizeX;
1047 const double yc = 0.;
1048 std::array<double, 3> centre = LocalToGlobal(xc, yc, zc, dcos, corner);
1049 // Assign element values.
1050 Element rect;
1051 rect.origin = centre;
1052 rect.lx = rectSizeX;
1053 rect.lz = elementSizeZ;
1054 rect.dA = rectSizeX * elementSizeZ;
1055 rect.dcos = dcos;
1056 // Boundary condition is applied at the origin of the rectangular
1057 // element coordinate system.
1058 rect.collocationPoint = centre;
1059 // Calculate the element vertices.
1060 const double hx = 0.5 * rectSizeX;
1061 const double hz = 0.5 * elementSizeZ;
1062 std::array<double, 3> p0 = LocalToGlobal(-hx, 0., -hz, dcos, centre);
1063 std::array<double, 3> p1 = LocalToGlobal( hx, 0., -hz, dcos, centre);
1064 std::array<double, 3> p2 = LocalToGlobal( hx, 0., hz, dcos, centre);
1065 std::array<double, 3> p3 = LocalToGlobal(-hx, 0., hz, dcos, centre);
1066 // Assign the vertices of the element.
1067 rect.xv = {p0[0], p1[0], p2[0], p3[0]};
1068 rect.yv = {p0[1], p1[1], p2[1], p3[1]};
1069 rect.zv = {p0[2], p1[2], p2[2], p3[2]};
1070 elements.push_back(std::move(rect));
1071 }
1072 }
1073 return true;
1074}
1075
1076bool ComponentNeBem3d::DiscretizeRectangle(const Primitive& primitive,
1077 const double targetSize, std::vector<Element>& elements) const {
1078
1079 std::array<double, 3> origin = {
1080 0.25 * std::accumulate(primitive.xv.begin(), primitive.xv.end(), 0.),
1081 0.25 * std::accumulate(primitive.yv.begin(), primitive.yv.end(), 0.),
1082 0.25 * std::accumulate(primitive.zv.begin(), primitive.zv.end(), 0.)};
1083
1084 // Determine the direction cosines.
1085 const double dx1 = primitive.xv[1] - primitive.xv[0];
1086 const double dy1 = primitive.yv[1] - primitive.yv[0];
1087 const double dz1 = primitive.zv[1] - primitive.zv[0];
1088 const double dx2 = primitive.xv[2] - primitive.xv[1];
1089 const double dy2 = primitive.yv[2] - primitive.yv[1];
1090 const double dz2 = primitive.zv[2] - primitive.zv[1];
1091 double lx = sqrt(dx1 * dx1 + dy1 * dy1 + dz1 * dz1);
1092 double lz = sqrt(dx2 * dx2 + dy2 * dy2 + dz2 * dz2);
1093 const std::array<double, 3> xu = {dx1 / lx, dy1 / lx, dz1 / lx};
1094 const std::array<double, 3> yu = {primitive.a, primitive.b, primitive.c};
1095 const std::array<double, 3> zu = CrossProduct(xu, yu);
1096 const std::array<std::array<double, 3>, 3> dcos = {xu, yu, zu};
1097
1098 // Determine the number of elements.
1099 unsigned int nx = NbOfSegments(lx, targetSize);
1100 unsigned int nz = NbOfSegments(lz, targetSize);
1101
1102 double elementSizeX = lx / nx;
1103 double elementSizeZ = lz / nz;
1104
1105 // Analyze the element aspect ratio.
1106 double ar = elementSizeX / elementSizeZ;
1107 if (ar > 10.) {
1108 // Reduce nz.
1109 elementSizeZ = 0.1 * elementSizeX;
1110 nz = std::max(static_cast<int>(lz / elementSizeZ), 1);
1111 elementSizeZ = lz / nz;
1112 }
1113 if (ar < 0.1) {
1114 // Reduce nx.
1115 elementSizeX = 0.1 * elementSizeZ;
1116 nx = std::max(static_cast<int>(lx / elementSizeX), 1);
1117 elementSizeX = lx / nx;
1118 }
1119
1120 const double dA = elementSizeX * elementSizeZ;
1121 for (unsigned int i = 0; i < nx; ++i) {
1122 // Centroid of the element.
1123 const double xav = -0.5 * lx + (i + 0.5) * elementSizeX;
1124 for (unsigned int k = 0; k < nz; ++k) {
1125 const double zav = -0.5 * lz + (k + 0.5) * elementSizeZ;
1126
1127 std::array<double, 3> centre = LocalToGlobal(xav, 0., zav, dcos, origin);
1128 Element element;
1129 element.origin = centre;
1130 element.lx = elementSizeX;
1131 element.lz = elementSizeZ;
1132 element.dA = dA;
1133 element.dcos = dcos;
1134 element.collocationPoint = centre;
1135 // Calculate the element vertices.
1136 const double hx = 0.5 * elementSizeX;
1137 const double hz = 0.5 * elementSizeZ;
1138 std::array<double, 3> p0 = LocalToGlobal(-hx, 0., -hz, dcos, centre);
1139 std::array<double, 3> p1 = LocalToGlobal( hx, 0., -hz, dcos, centre);
1140 std::array<double, 3> p2 = LocalToGlobal( hx, 0., hz, dcos, centre);
1141 std::array<double, 3> p3 = LocalToGlobal(-hx, 0., hz, dcos, centre);
1142 // Assign the vertices of the element.
1143 element.xv = {p0[0], p1[0], p2[0], p3[0]};
1144 element.yv = {p0[1], p1[1], p2[1], p3[1]};
1145 element.zv = {p0[2], p1[2], p2[2], p3[2]};
1146 elements.push_back(std::move(element));
1147 }
1148 }
1149 return true;
1150}
1151
1152unsigned int ComponentNeBem3d::NbOfSegments(const double length,
1153 const double target) const {
1154
1155 // Check whether the length of the primitive is long enough.
1156 if (length < MinDist) return 1;
1157 unsigned int n = static_cast<unsigned int>(length / target);
1158 if (n < m_minNbElementsOnLength) {
1159 // Need to have a minimum number of elements per primitive...
1160 n = m_minNbElementsOnLength;
1161 if (length < n * MinDist) {
1162 // ...which may not be possible if the length is small.
1163 n = static_cast<unsigned int>(length / MinDist);
1164 if (n < 1) {
1165 // However, it is necessary to have at least one element!
1166 n = 1;
1167 }
1168 }
1169 }
1170 return std::min(n, m_maxNbElementsOnLength);
1171}
1172
1173bool ComponentNeBem3d::EliminateOverlaps(const Panel& panel1,
1174 const Panel& panel2,
1175 std::vector<Panel>& panelsOut,
1176 std::vector<int>& itypo) {
1177 // *-----------------------------------------------------------------------
1178 // * PLAOVL - Isolates the parts of plane 1 that are not hidden by 2.
1179 // *-----------------------------------------------------------------------
1180
1181 const auto& xp1 = panel1.xv;
1182 const auto& yp1 = panel1.yv;
1183 const auto& zp1 = panel1.zv;
1184 const auto& xp2 = panel2.xv;
1185 const auto& yp2 = panel2.yv;
1186 const auto& zp2 = panel2.zv;
1187 // If the size of either is less than 3, simply return.
1188 if (xp1.size() <= 2 || xp2.size() <= 2) {
1189 return true;
1190 }
1191 // Compute the various tolerances.
1192 const double xmin1 = *std::min_element(std::begin(xp1), std::end(xp1));
1193 const double ymin1 = *std::min_element(std::begin(yp1), std::end(yp1));
1194 const double xmax1 = *std::max_element(std::begin(xp1), std::end(xp1));
1195 const double ymax1 = *std::max_element(std::begin(yp1), std::end(yp1));
1196
1197 const double xmin2 = *std::min_element(std::begin(xp2), std::end(xp2));
1198 const double ymin2 = *std::min_element(std::begin(yp2), std::end(yp2));
1199 const double xmax2 = *std::max_element(std::begin(xp2), std::end(xp2));
1200 const double ymax2 = *std::max_element(std::begin(yp2), std::end(yp2));
1201
1202 const double xmin = std::min(xmin1, xmin2);
1203 const double ymin = std::min(ymin1, ymin2);
1204 const double xmax = std::max(xmax1, xmax2);
1205 const double ymax = std::max(ymax1, ymax2);
1206
1207 const double epsx = 1.e-6 * std::max(std::abs(xmax), std::abs(xmin));
1208 const double epsy = 1.e-6 * std::max(std::abs(ymax), std::abs(ymin));
1209
1210 const double zsum1 = std::accumulate(std::begin(zp1), std::end(zp1), 0.);
1211 const double zsum2 = std::accumulate(std::begin(zp2), std::end(zp2), 0.);
1212 const double zmean = (zsum1 + zsum2) / (zp1.size() + zp2.size());
1213
1214 std::array<std::vector<double>, 2> xl;
1215 std::array<std::vector<double>, 2> yl;
1216 std::array<std::vector<int>, 2> flags;
1217 std::array<std::vector<double>, 2> qs;
1218 // Establish the list of special points around polygon 1.
1219 AddPoints(xp1, yp1, xp2, yp2, xl[0], yl[0], flags[0], qs[0], epsx, epsy);
1220 // Establish the list of special points around polygon 2.
1221 AddPoints(xp2, yp2, xp1, yp1, xl[1], yl[1], flags[1], qs[1], epsx, epsy);
1222
1223 bool ok = true;
1224 // Look up the cross-links: from plane 1 (2) to plane 2 (1).
1225 std::array<std::vector<int>, 2> links;
1226 for (unsigned int ic = 0; ic < 2; ++ic) {
1227 const unsigned int n1 = xl[ic].size();
1228 links[ic].assign(n1, -1);
1229 const unsigned int jc = ic == 0 ? 1 : 0;
1230 const unsigned int n2 = xl[jc].size();
1231 for (unsigned int i = 0; i < n1; ++i) {
1232 unsigned int nFound = 0;
1233 for (unsigned int j = 0; j < n2; ++j) {
1234 if (fabs(xl[ic][i] - xl[jc][j]) < epsx &&
1235 fabs(yl[ic][i] - yl[jc][j]) < epsy) {
1236 ++nFound;
1237 links[ic][i] = j;
1238 }
1239 }
1240 if (nFound == 0 && (flags[ic][i] == 2 || flags[ic][i] == 3)) {
1241 std::cerr << m_className << "::EliminateOverlaps: "
1242 << "Warning. Expected match not found (" << ic << "-"
1243 << jc << ").\n";
1244 links[ic][i] = -1;
1245 ok = false;
1246 } else if (nFound > 1) {
1247 std::cerr << m_className << "::EliminateOverlaps: "
1248 << "Warning. More than 1 match found (" << ic << "-"
1249 << jc << ").\n";
1250 links[ic][i] = -1;
1251 ok = false;
1252 }
1253 }
1254 }
1255
1256 // List the points for debugging.
1257 if (m_debug) {
1258 for (unsigned int j = 0; j < 2; ++j) {
1259 std::cout << " Polygon " << j << "\n "
1260 << " No Type x y Q links\n";
1261 const unsigned int n = xl[j].size();
1262 for (unsigned int i = 0; i < n; ++i) {
1263 printf(" %3d %5d %13.6f %13.6f %5.3f %3d\n", i, flags[j][i],
1264 xl[j][i], yl[j][i], qs[j][i], links[j][i]);
1265 }
1266 }
1267 }
1268 if (!ok) return false;
1269
1270 for (unsigned int ic = 0; ic < 2; ++ic) {
1271 // See whether all of 1 (2) is inside 2 (1).
1272 bool allInside = true;
1273 const unsigned int np = xl[ic].size();
1274 for (unsigned int i = 0; i < np; ++i) {
1275 if (flags[ic][i] != 1) {
1276 allInside = false;
1277 break;
1278 }
1279 bool inside = false, edge = false;
1280 if (ic == 0) {
1281 Polygon::Inside(xp2, yp2, xl[ic][i], yl[ic][i], inside, edge);
1282 } else {
1283 Polygon::Inside(xp1, yp1, xl[ic][i], yl[ic][i], inside, edge);
1284 }
1285 if (!(inside || edge)) {
1286 allInside = false;
1287 break;
1288 }
1289 }
1290 if (allInside) {
1291 // Apparently 1 (2) really is fully inside 2 (1).
1292 if (ic == 0) {
1293 if (m_debug) std::cout << " Curve 0 fully inside 1.\n";
1294 // Write out curve 1.
1295 panelsOut.push_back(panel1);
1296 } else {
1297 if (m_debug) std::cout << " Curve 1 fully inside 0.\n";
1298 // Write out curve 2.
1299 panelsOut.push_back(panel2);
1300 }
1301 panelsOut.back().zv.assign(panelsOut.back().xv.size(), zmean);
1302 itypo.push_back(3);
1303 std::vector<Panel> newPanels;
1304 if (ic == 0) {
1305 if (!TraceEnclosed(xl[0], yl[0], xl[1], yl[1], panel2, newPanels)) {
1306 return false;
1307 }
1308 } else {
1309 if (!TraceEnclosed(xl[1], yl[1], xl[0], yl[0], panel1, newPanels)) {
1310 return false;
1311 }
1312 }
1313 for (auto& panel : newPanels) {
1314 panel.zv.assign(panel.xv.size(), zmean);
1315 if (ic == 0) {
1316 itypo.push_back(2);
1317 } else {
1318 itypo.push_back(1);
1319 }
1320 }
1321 panelsOut.insert(panelsOut.end(), newPanels.begin(), newPanels.end());
1322 return true;
1323 }
1324 }
1325
1326 for (unsigned int ic = 0; ic < 2; ++ic) {
1327 std::vector<Panel> newPanels;
1328 const unsigned int n = xl[ic].size();
1329 // Identify the parts of 1 (2) that are not overlapped, first mark.
1330 std::vector<bool> mark(n, false);
1331 bool done = false;
1332 while (!done) {
1333 if (m_debug) {
1334 std::cout << " Searching for starting point on " << ic << ".\n";
1335 }
1336 done = true;
1337 // Try and find a new starting point
1338 for (unsigned int i = 0; i < n; ++i) {
1339 const unsigned int ii = NextPoint(i, n);
1340 // Skip parts already processed.
1341 if (mark[i] || mark[ii]) continue;
1342 // Skip if mid point is inside other volume.
1343 bool inside = false, edge = false;
1344 const double xm = 0.5 * (xl[ic][i] + xl[ic][ii]);
1345 const double ym = 0.5 * (yl[ic][i] + yl[ic][ii]);
1346 if (ic == 0) {
1347 Polygon::Inside(xp2, yp2, xm, ym, inside, edge);
1348 } else {
1349 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
1350 }
1351 if (inside || edge) continue;
1352 // Found one.
1353 done = false;
1354
1355 if (ic == 0) {
1356 // Trace this part of 1 outside 2.
1357 TraceNonOverlap(xp1, yp1, xl[0], yl[0], xl[1], yl[1], flags[0],
1358 flags[1], links[0], links[1], mark, i, panel1,
1359 newPanels);
1360 } else {
1361 // Trace this part of 2 outside 1.
1362 TraceNonOverlap(xp2, yp2, xl[1], yl[1], xl[0], yl[0], flags[1],
1363 flags[0], links[1], links[0], mark, i, panel2,
1364 newPanels);
1365 }
1366 break;
1367 }
1368 }
1369 for (auto& panel : newPanels) {
1370 panel.zv.assign(panel.xv.size(), zmean);
1371 itypo.push_back(ic + 1);
1372 }
1373 panelsOut.insert(panelsOut.end(), newPanels.begin(), newPanels.end());
1374 if (m_debug) {
1375 std::cout << " No further non-overlapped areas of " << ic << ".\n";
1376 }
1377 }
1378
1379 // Look for the overlapped parts.
1380 std::vector<Panel> newPanels;
1381 const unsigned int n1 = xl[0].size();
1382 std::vector<bool> mark1(n1, false);
1383 bool done = false;
1384 while (!done) {
1385 done = true;
1386 if (m_debug) {
1387 std::cout << " Searching for starting point on overlap.\n";
1388 }
1389 for (unsigned int i = 0; i < n1; ++i) {
1390 // Skip points already processed.
1391 if (mark1[i]) continue;
1392 // Skip if not an edge point on both 1 and 2 or internal in 2.
1393 int ip1 = i;
1394 int ip2 = links[0][ip1];
1395 if (ip2 < 0 || flags[0][ip1] == 1) {
1396 bool inside = false, edge = false;
1397 Polygon::Inside(xp2, yp2, xl[0][ip1], yl[0][ip1], inside, edge);
1398 if (!(inside || edge)) continue;
1399 } else if (flags[1][ip2] == 1) {
1400 continue;
1401 }
1402 // Found one.
1403 done = false;
1404 TraceOverlap(xp1, yp1, xp2, yp2, xl[0], yl[0], xl[1], yl[1], flags[0],
1405 links[0], links[1], mark1, ip1, ip2, panel1, newPanels);
1406 break;
1407 }
1408 }
1409 for (auto& panel : newPanels) {
1410 panel.zv.assign(panel.xv.size(), zmean);
1411 itypo.push_back(3);
1412 }
1413 panelsOut.insert(panelsOut.end(), newPanels.begin(), newPanels.end());
1414 // Finished
1415 if (m_debug) std::cout << " No further overlapped areas.\n";
1416 return true;
1417}
1418
1419bool ComponentNeBem3d::TraceEnclosed(const std::vector<double>& xl1,
1420 const std::vector<double>& yl1,
1421 const std::vector<double>& xl2,
1422 const std::vector<double>& yl2,
1423 const Panel& panel2,
1424 std::vector<Panel>& panelsOut) const {
1425 const int n1 = xl1.size();
1426 const int n2 = xl2.size();
1427 // Find 2 non-crossing connections: JP1-JP2 and KP1-KP2.
1428 unsigned int nFound = 0;
1429 int jp1 = 0, jp2 = 0;
1430 int kp1 = 0, kp2 = 0;
1431 for (int ip1 = 0; ip1 < n1; ++ip1) {
1432 const double x1 = xl1[ip1];
1433 const double y1 = yl1[ip1];
1434 for (int ip2 = 0; ip2 < n2; ++ip2) {
1435 if (nFound > 0 && ip2 == jp2) continue;
1436 const double x2 = xl2[ip2];
1437 const double y2 = yl2[ip2];
1438 bool cross = false;
1439 for (int k = 0; k < n1; ++k) {
1440 const int kk = NextPoint(k, n1);
1441 if (k == ip1 || kk == ip1) continue;
1442 double xc = 0., yc = 0.;
1443 cross =
1444 Crossing(x1, y1, x2, y2, xl1[k], yl1[k], xl1[kk], yl1[kk], xc, yc);
1445 if (cross) break;
1446 }
1447 if (cross) continue;
1448 if (m_debug) std::cout << " No crossing with 1.\n";
1449 for (int k = 0; k < n2; ++k) {
1450 const int kk = NextPoint(k, n2);
1451 if (k == ip2 || kk == ip2) continue;
1452 double xc = 0., yc = 0.;
1453 cross =
1454 Crossing(x1, y1, x2, y2, xl2[k], yl2[k], xl2[kk], yl2[kk], xc, yc);
1455 if (cross) break;
1456 }
1457 if (cross) continue;
1458 if (nFound == 0) {
1459 jp1 = ip1;
1460 jp2 = ip2;
1461 if (m_debug) {
1462 std::cout << " 1st junction: " << jp1 << ", " << jp2 << ".\n";
1463 }
1464 ++nFound;
1465 break;
1466 } else {
1467 kp1 = ip1;
1468 kp2 = ip2;
1469 double xc = 0., yc = 0.;
1470 cross = Crossing(x1, y1, x2, y2, xl1[jp1], yl1[jp1], xl2[jp2], yl2[jp2],
1471 xc, yc);
1472 if (!cross) {
1473 if (m_debug) {
1474 std::cout << " 2nd junction: " << kp1 << ", " << kp2 << ".\n";
1475 }
1476 ++nFound;
1477 break;
1478 }
1479 }
1480 }
1481 if (nFound > 1) break;
1482 }
1483 if (nFound < 2) {
1484 std::cerr << m_className << "::TraceEnclosed: Found no cut-out.\n";
1485 return false;
1486 }
1487
1488 // Create part 1 of area 2.
1489 std::vector<double> xpl;
1490 std::vector<double> ypl;
1491 if (m_debug) std::cout << " Creating part 1 of area 2.\n";
1492 for (int ip1 = jp1; ip1 <= kp1; ++ip1) {
1493 if (m_debug) std::cout << " Adding " << ip1 << " on 1.\n";
1494 xpl.push_back(xl1[ip1]);
1495 ypl.push_back(yl1[ip1]);
1496 }
1497 // Try one way.
1498 int imax = jp2 < kp2 ? jp2 + n2 : jp2;
1499 int dir = +1;
1500 for (int i = kp2; i <= imax; ++i) {
1501 int ip2 = i % n2;
1502 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
1503 xpl.push_back(xl2[ip2]);
1504 ypl.push_back(yl2[ip2]);
1505 }
1506 // Check for undesirable crossings.
1507 bool ok = true;
1508 for (int ip1 = 0; ip1 < n1; ++ip1) {
1509 if (ip1 == jp1 || ip1 == kp1) continue;
1510 bool inside = false, edge = false;
1511 Polygon::Inside(xpl, ypl, xl1[ip1], yl1[ip1], inside, edge);
1512 if (inside) {
1513 ok = false;
1514 break;
1515 }
1516 }
1517 if (!ok) {
1518 // Use the other way if this failed
1519 if (m_debug) std::cout << " Trying the other direction.\n";
1520 xpl.resize(kp1 - jp1 + 1);
1521 ypl.resize(kp1 - jp1 + 1);
1522 imax = jp2 < kp2 ? kp2 : kp2 + n2;
1523 dir = -1;
1524 for (int i = imax; i >= jp2; --i) {
1525 const int ip2 = i % n2;
1526 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
1527 xpl.push_back(xl2[ip2]);
1528 ypl.push_back(yl2[ip2]);
1529 }
1530 }
1531
1532 // Save this part.
1533 Panel newPanel1 = panel2;
1534 newPanel1.xv = xpl;
1535 newPanel1.yv = ypl;
1536 panelsOut.push_back(std::move(newPanel1));
1537 if (m_debug) std::cout << " Part 1 has " << xpl.size() << " nodes.\n";
1538
1539 // Create part 2 of area 2.
1540 xpl.clear();
1541 ypl.clear();
1542 if (m_debug) std::cout << " Creating part 2 of area 2.\n";
1543 imax = jp1 + n1;
1544 for (int i = kp1; i <= imax; ++i) {
1545 const int ip1 = i % n1;
1546 if (m_debug) std::cout << " Adding " << ip1 << " on 1.\n";
1547 xpl.push_back(xl1[ip1]);
1548 ypl.push_back(yl1[ip1]);
1549 }
1550 // Add the part over area 2.
1551 if (dir == -1) {
1552 imax = jp2 > kp2 ? jp2 : jp2 + n2;
1553 for (int i = imax; i >= kp2; --i) {
1554 const int ip2 = i % n2;
1555 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
1556 xpl.push_back(xl2[ip2]);
1557 ypl.push_back(yl2[ip2]);
1558 }
1559 } else {
1560 imax = jp2 > kp2 ? kp2 + n2 : kp2;
1561 for (int i = jp2; i <= imax; ++i) {
1562 const int ip2 = i % n2;
1563 if (m_debug) std::cout << " Adding " << ip2 << " on 2.\n";
1564 xpl.push_back(xl2[ip2]);
1565 ypl.push_back(yl2[ip2]);
1566 }
1567 }
1568 // Save this part.
1569 Panel newPanel2 = panel2;
1570 newPanel2.xv = xpl;
1571 newPanel2.yv = ypl;
1572 panelsOut.push_back(std::move(newPanel2));
1573 if (m_debug) std::cout << " Part 1 has " << xpl.size() << " nodes.\n";
1574 return true;
1575}
1576
1577void ComponentNeBem3d::TraceNonOverlap(
1578 const std::vector<double>& xp1, const std::vector<double>& yp1,
1579 const std::vector<double>& xl1, const std::vector<double>& yl1,
1580 const std::vector<double>& xl2, const std::vector<double>& yl2,
1581 const std::vector<int>& flags1, const std::vector<int>& flags2,
1582 const std::vector<int>& links1, const std::vector<int>& links2,
1583 std::vector<bool>& mark1, int ip1, const Panel& panel1,
1584 std::vector<Panel>& panelsOut) const {
1585 const unsigned int n1 = xl1.size();
1586 const unsigned int n2 = xl2.size();
1587
1588 // Remember the starting point.
1589 const int is1 = ip1;
1590 std::vector<double> xpl;
1591 std::vector<double> ypl;
1592 // Add the starting point.
1593 xpl.push_back(xl1[ip1]);
1594 ypl.push_back(yl1[ip1]);
1595 mark1[ip1] = true;
1596 if (m_debug) {
1597 std::cout << " Start from point " << ip1 << " on curve 1.\n";
1598 }
1599 // Next point.
1600 ip1 = NextPoint(ip1, n1);
1601 xpl.push_back(xl1[ip1]);
1602 ypl.push_back(yl1[ip1]);
1603 mark1[ip1] = true;
1604 if (m_debug) {
1605 std::cout << " Next point is " << ip1 << " on curve 1.\n";
1606 }
1607
1608 // Keep track of the curve we are currently following.
1609 unsigned int il = 1;
1610 // Direction flag (-1: backward, 0: not set, 1: forward).
1611 int dir = 0;
1612 // End-of-curve flag.
1613 bool eoc = false;
1614 int ip2 = 0;
1615 while (!eoc) {
1616 if (il == 1 && flags1[std::max(ip1, 0)] == 1) {
1617 // On curve 1 and not on the edge of curve 2?
1618 ip1 = NextPoint(ip1, n1);
1619 if (ip1 == is1) {
1620 eoc = true;
1621 continue;
1622 }
1623 mark1[ip1] = true;
1624 xpl.push_back(xl1[ip1]);
1625 ypl.push_back(yl1[ip1]);
1626 if (m_debug) {
1627 std::cout << " Went to point " << ip1 << " on curve 1.\n";
1628 }
1629 } else if (il == 1) {
1630 // On curve 1 and on the edge of curve 2?
1631 ip2 = links1[ip1];
1632 bool added = false;
1633 if (dir == +1 || dir == 0) {
1634 const double xm = 0.5 * (xl2[ip2] + xl2[NextPoint(ip2, n2)]);
1635 const double ym = 0.5 * (yl2[ip2] + yl2[NextPoint(ip2, n2)]);
1636 bool inside = false, edge = false;
1637 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
1638 if (inside) {
1639 ip2 = NextPoint(ip2, n2);
1640 il = 2;
1641 dir = +1;
1642 ip1 = links2[ip2];
1643 if (ip1 == is1) {
1644 eoc = true;
1645 continue;
1646 } else if (ip1 >= 0) {
1647 mark1[ip1] = true;
1648 }
1649 xpl.push_back(xl2[ip2]);
1650 ypl.push_back(yl2[ip2]);
1651 added = true;
1652 if (m_debug) {
1653 std::cout << " Added point " << ip2 << " along 2 +.\n";
1654 }
1655 }
1656 }
1657 if (dir == -1 || dir == 0) {
1658 const double xm = 0.5 * (xl2[ip2] + xl2[PrevPoint(ip2, n2)]);
1659 const double ym = 0.5 * (yl2[ip2] + yl2[PrevPoint(ip2, n2)]);
1660 bool inside = false, edge = false;
1661 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
1662 if (inside) {
1663 ip2 = PrevPoint(ip2, n2);
1664 il = 2;
1665 dir = -1;
1666 ip1 = links2[ip2];
1667 if (ip1 == is1) {
1668 eoc = true;
1669 continue;
1670 } else if (ip1 >= 0) {
1671 mark1[ip1] = true;
1672 }
1673 xpl.push_back(xl2[ip2]);
1674 ypl.push_back(yl2[ip2]);
1675 added = true;
1676 if (m_debug) {
1677 std::cout << " Added point " << ip2 << " along 2 -\n";
1678 }
1679 }
1680 }
1681 if (!added) {
1682 ip1 = NextPoint(ip1, n1);
1683 if (ip1 == is1) {
1684 eoc = true;
1685 continue;
1686 } else if (ip1 >= 0) {
1687 mark1[ip1] = true;
1688 }
1689 xpl.push_back(xl1[ip1]);
1690 ypl.push_back(yl1[ip1]);
1691 if (m_debug) std::cout << " Continued over 1.\n";
1692 }
1693 } else if (il == 2 && flags2[std::max(ip2, 0)] == 1) {
1694 // On curve 2 normal vertex (outside 1 hopefully).
1695 ip2 = dir > 0 ? NextPoint(ip2, n2) : PrevPoint(ip2, n2);
1696 ip1 = links2[ip2];
1697 if (ip1 == is1) {
1698 eoc = true;
1699 continue;
1700 } else if (ip1 >= 0) {
1701 mark1[ip1] = true;
1702 }
1703 xpl.push_back(xl2[ip2]);
1704 ypl.push_back(yl2[ip2]);
1705 if (m_debug) {
1706 std::cout << " Went to point " << ip2 << " on 2.\n";
1707 }
1708 } else if (il == 2) {
1709 // On curve 2 and on edge of 1.
1710 ip1 = links2[ip2];
1711 ip1 = NextPoint(ip1, n1);
1712 il = 1;
1713 if (ip1 == is1) {
1714 eoc = true;
1715 continue;
1716 }
1717 xpl.push_back(xl1[ip1]);
1718 ypl.push_back(yl1[ip1]);
1719 if (m_debug) std::cout << " Resumed 1 at point " << ip1 << ".\n";
1720 } else {
1721 // Other cases should not occur.
1722 std::cerr << m_className << "::TraceNonOverlap: Unexpected case.\n";
1723 return;
1724 }
1725 }
1726
1727 Panel newPanel = panel1;
1728 newPanel.xv = xpl;
1729 newPanel.yv = ypl;
1730 panelsOut.push_back(std::move(newPanel));
1731 if (m_debug) {
1732 std::cout << " End of curve reached, " << xpl.size() << " points.\n";
1733 }
1734}
1735
1736void ComponentNeBem3d::TraceOverlap(
1737 const std::vector<double>& xp1, const std::vector<double>& yp1,
1738 const std::vector<double>& xp2, const std::vector<double>& yp2,
1739 const std::vector<double>& xl1, const std::vector<double>& yl1,
1740 const std::vector<double>& xl2, const std::vector<double>& yl2,
1741 const std::vector<int>& flags1, const std::vector<int>& links1,
1742 const std::vector<int>& links2, std::vector<bool>& mark1, int ip1, int ip2,
1743 const Panel& panel1, std::vector<Panel>& panelsOut) const {
1744 int ip1L = -1;
1745 int ip1LL = -1;
1746
1747 const unsigned int n1 = xl1.size();
1748 const unsigned int n2 = xl2.size();
1749
1750 // Remember the starting points.
1751 const int is1 = ip1;
1752 const int is2 = ip2;
1753 std::vector<double> xpl;
1754 std::vector<double> ypl;
1755 xpl.push_back(xl1[ip1]);
1756 ypl.push_back(yl1[ip1]);
1757 mark1[ip1] = true;
1758
1759 // Keep track of the curve we are currently following.
1760 unsigned int il = 1;
1761 // Direction flag (-1: backward, 0: not set, 1: forward).
1762 int dir = 0;
1763 // End-of-curve flag.
1764 bool eoc = false;
1765
1766 if (m_debug) {
1767 std::cout << " Start from point " << ip1 << " on curve " << il << "\n";
1768 }
1769 while (!eoc) {
1770 ip1LL = ip1L;
1771 ip1L = ip1;
1772 if (il == 1) {
1773 // We are on curve 1. Move to the next point.
1774 const int ii = NextPoint(ip1, n1);
1775 // Maybe finished over line 1?
1776 if (ii == is1) {
1777 eoc = true;
1778 continue;
1779 }
1780 // See whether the next point of 1 is on the edge or inside of 2.
1781 bool inside = false, edge = false;
1782 if (links1[ii] >= 0) {
1783 edge = true;
1784 } else if (flags1[ii] == 1) {
1785 Polygon::Inside(xp2, yp2, xl1[ii], yl1[ii], inside, edge);
1786 }
1787 // If it is, check that it doesn't leave 2 at any stage.
1788 if (inside || edge) {
1789 const double xm = 0.5 * (xl1[ip1] + xl1[ii]);
1790 const double ym = 0.5 * (yl1[ip1] + yl1[ii]);
1791 Polygon::Inside(xp2, yp2, xm, ym, inside, edge);
1792 }
1793 // If it is, continue over 1.
1794 if (inside || edge) {
1795 ip1 = ii;
1796 if (m_debug) {
1797 std::cout << " Continued to point " << ip1 << " on "
1798 << il << "\n";
1799 }
1800 xpl.push_back(xl1[ip1]);
1801 ypl.push_back(yl1[ip1]);
1802 mark1[ip1] = true;
1803 continue;
1804 }
1805 // Else we have to continue over 2, ensure we really are on curve 2.
1806 ip2 = links1[ip1];
1807 if (ip2 < 0) {
1808 std::cerr << m_className << "::TraceOverlap: "
1809 << "No point 2 reference found; abandoned.\n";
1810 return;
1811 }
1812 // Impose a direction on 2 to avoid returning.
1813 if (dir == 0) {
1814 if (links2[NextPoint(ip2, n2)] == ip1LL &&
1815 links2[PrevPoint(ip2, n2)] == ip1LL) {
1816 std::cerr << m_className << "::TraceOverlap: "
1817 << "Both 2+ and 2- return on 1; not stored.\n";
1818 return;
1819 } else if (links2[NextPoint(ip2, n2)] == ip1LL) {
1820 if (m_debug) {
1821 std::cout << " 2+ is a return to previous point on 1.\n";
1822 }
1823 dir = -1;
1824 } else if (links2[PrevPoint(ip2, n2)] == ip1LL) {
1825 if (m_debug) {
1826 std::cout << " 2- is a return to previous point on 1.\n";
1827 }
1828 dir = +1;
1829 } else {
1830 if (m_debug) std::cout << " Both ways are OK.\n";
1831 }
1832 }
1833 // If not, try to continue over 2 in the + direction.
1834 if (dir == +1 || dir == 0) {
1835 ip2 = NextPoint(ip2, n2);
1836 if (ip2 == is2) {
1837 if (m_debug) std::cout << " Return to start over 2+.\n";
1838 eoc = true;
1839 continue;
1840 }
1841 Polygon::Inside(xp1, yp1, xl2[ip2], yl2[ip2], inside, edge);
1842 if (inside || edge) {
1843 if (m_debug) {
1844 std::cout << " Going to 2+ (point " << ip2 << " of 2).\n";
1845 }
1846 xpl.push_back(xl2[ip2]);
1847 ypl.push_back(yl2[ip2]);
1848 dir = +1;
1849 if (links2[ip2] >= 0) {
1850 ip1 = links2[ip2];
1851 mark1[ip1] = true;
1852 il = 1;
1853 if (m_debug) {
1854 std::cout << " This point is also on curve 1: "
1855 << ip1 << "\n";
1856 }
1857 } else {
1858 il = 2;
1859 }
1860 continue;
1861 }
1862 // Continuing in the + direction didn't work so go back a step.
1863 ip2 = PrevPoint(ip2, n2);
1864 }
1865 // Or if this still fails, try 2 in the - direction.
1866 if (dir == -1 || dir == 0) {
1867 ip2 = PrevPoint(ip2, n2);
1868 if (ip2 == is2) {
1869 if (m_debug) std::cout << " Return to start over 2-\n";
1870 eoc = true;
1871 continue;
1872 }
1873 Polygon::Inside(xp1, yp1, xl2[ip2], yl2[ip2], inside, edge);
1874 if (inside || edge) {
1875 if (m_debug) {
1876 std::cout << " Going to 2- (point " << ip2 << " of 2).\n";
1877 }
1878 xpl.push_back(xl2[ip2]);
1879 ypl.push_back(yl2[ip2]);
1880 dir = -1;
1881 if (links2[ip2] >= 0) {
1882 ip1 = links2[ip2];
1883 mark1[ip1] = true;
1884 il = 1;
1885 if (m_debug) {
1886 std::cout << " This point is also on 1: " << ip1 << ".\n";
1887 }
1888 } else {
1889 il = 2;
1890 }
1891 continue;
1892 }
1893 }
1894 // Should not get here.
1895 if (m_debug) std::cout << " Dead end.\n";
1896 return;
1897 } else if (il == 2) {
1898 // We are on curve 2. Ensure the direction is set.
1899 if (dir == 0) {
1900 std::cerr << m_className << "::TraceOverlap: "
1901 << "Direction not set; abandoned.\n";
1902 return;
1903 }
1904 // Move to the next point.
1905 ip2 = dir > 0 ? NextPoint(ip2, n2) : PrevPoint(ip2, n2);
1906 // Maybe finished over line 2?
1907 if (ip2 == is2) {
1908 // Reached the end.
1909 eoc = true;
1910 continue;
1911 }
1912 // Next step over 2.
1913 if (m_debug) {
1914 std::cout << " Stepped over 2 to point " << ip2 << " of 2.\n";
1915 }
1916 xpl.push_back(xl2[ip2]);
1917 ypl.push_back(yl2[ip2]);
1918 if (links2[ip2] >= 0) {
1919 ip1 = links2[ip2];
1920 mark1[ip1] = true;
1921 il = 1;
1922 if (m_debug) {
1923 std::cout << " This point is also on curve 1: " << ip1 << ".\n";
1924 }
1925 } else {
1926 il = 2;
1927 }
1928 }
1929 }
1930
1931 if (xpl.size() <= 2) {
1932 if (m_debug) std::cout << " Too few points.\n";
1933 } else {
1934 Panel newPanel = panel1;
1935 newPanel.xv = xpl;
1936 newPanel.yv = ypl;
1937 panelsOut.push_back(std::move(newPanel));
1938 }
1939
1940 if (m_debug) {
1941 std::cout << " End of curve reached, " << xpl.size() << " points.\n";
1942 }
1943}
1944
1945bool ComponentNeBem3d::MakePrimitives(const Panel& panelIn,
1946 std::vector<Panel>& panelsOut) const {
1947 // *-----------------------------------------------------------------------
1948 // * PLATRC - Cuts a polygon into right-angled triangles.
1949 // *-----------------------------------------------------------------------
1950
1951 // Establish tolerances.
1952 // TODO! Class member?
1953 const double epsang = 1.e-6; // BEMEPA
1954
1955 if (panelIn.xv.empty() || panelIn.yv.empty() || panelIn.zv.empty()) {
1956 return false;
1957 }
1958 // Determine the mean z value.
1959 const double zsum =
1960 std::accumulate(std::begin(panelIn.zv), std::end(panelIn.zv), 0.);
1961 const double zmean = zsum / panelIn.zv.size();
1962
1963 std::vector<Panel> stack;
1964 stack.push_back(panelIn);
1965 stack.back().zv.clear();
1966 for (unsigned int k = 0; k < stack.size(); ++k) {
1967 // Next polygon.
1968 const auto& xp1 = stack[k].xv;
1969 const auto& yp1 = stack[k].yv;
1970 const unsigned int np = xp1.size();
1971 if (m_debug) {
1972 std::cout << " Polygon " << k << " with " << np << " nodes.\n";
1973 }
1974 if (np <= 2) {
1975 // Too few nodes.
1976 if (m_debug) std::cout << " Too few points.\n";
1977 continue;
1978 }
1979
1980 // See whether this is a right-angled triangle.
1981 if (np == 3) {
1982 const double x12 = xp1[0] - xp1[1];
1983 const double y12 = yp1[0] - yp1[1];
1984 const double x23 = xp1[1] - xp1[2];
1985 const double y23 = yp1[1] - yp1[2];
1986 const double x31 = xp1[2] - xp1[0];
1987 const double y31 = yp1[2] - yp1[0];
1988 const double x32 = xp1[2] - xp1[1];
1989 const double y32 = yp1[2] - yp1[1];
1990 const double x13 = xp1[0] - xp1[2];
1991 const double y13 = yp1[0] - yp1[2];
1992 const double x21 = xp1[1] - xp1[0];
1993 const double y21 = yp1[1] - yp1[0];
1994 const double s12 = x12 * x12 + y12 * y12;
1995 const double s32 = x32 * x32 + y32 * y32;
1996 const double s13 = x13 * x13 + y13 * y13;
1997 const double s23 = x23 * x23 + y23 * y23;
1998 const double s31 = x31 * x31 + y31 * y31;
1999 const double s21 = x21 * x21 + y21 * y21;
2000 if (fabs(x12 * x32 + y12 * y32) < epsang * sqrt(s12 * s32)) {
2001 if (m_debug) {
2002 std::cout << " Right-angled triangle node 2 - done.\n";
2003 }
2004 panelsOut.push_back(stack[k]);
2005 continue;
2006 } else if (fabs(x13 * x23 + y13 * y23) < epsang * sqrt(s13 * s23)) {
2007 if (m_debug) {
2008 std::cout << " Right-angled triangle node 3 - rearrange.\n";
2009 }
2010 Panel panel = stack[k];
2011 panel.xv = {xp1[1], xp1[2], xp1[0]};
2012 panel.yv = {yp1[1], yp1[2], yp1[0]};
2013 panelsOut.push_back(std::move(panel));
2014 continue;
2015 } else if (fabs(x31 * x21 + y31 * y21) < epsang * sqrt(s31 * s21)) {
2016 if (m_debug) {
2017 std::cout << " Right-angled triangle node 1 - rearrange.\n";
2018 }
2019 Panel panel = stack[k];
2020 panel.xv = {xp1[2], xp1[0], xp1[1]};
2021 panel.yv = {yp1[2], yp1[0], yp1[1]};
2022 panelsOut.push_back(std::move(panel));
2023 continue;
2024 }
2025 }
2026 // See whether this is a rectangle.
2027 if (np == 4) {
2028 const double x12 = xp1[0] - xp1[1];
2029 const double y12 = yp1[0] - yp1[1];
2030 const double x23 = xp1[1] - xp1[2];
2031 const double y23 = yp1[1] - yp1[2];
2032 const double x34 = xp1[2] - xp1[3];
2033 const double y34 = yp1[2] - yp1[3];
2034 const double x43 = xp1[3] - xp1[2];
2035 const double y43 = yp1[3] - yp1[2];
2036 const double x14 = xp1[0] - xp1[3];
2037 const double y14 = yp1[0] - yp1[3];
2038 const double x32 = xp1[2] - xp1[1];
2039 const double y32 = yp1[2] - yp1[1];
2040
2041 const double s12 = x12 * x12 + y12 * y12;
2042 const double s23 = x23 * x23 + y23 * y23;
2043 const double s32 = x32 * x32 + y32 * y32;
2044 const double s34 = x34 * x34 + y34 * y34;
2045 const double s43 = x43 * x43 + y43 * y43;
2046 const double s14 = x14 * x14 + y14 * y14;
2047 if (fabs(x12 * x32 + y12 * y32) < epsang * sqrt(s12 * s32) &&
2048 fabs(x23 * x43 + y23 * y43) < epsang * sqrt(s23 * s43) &&
2049 fabs(x14 * x34 + y14 * y34) < epsang * sqrt(s14 * s34)) {
2050 if (m_debug) std::cout << " Rectangle.\n";
2051 panelsOut.push_back(stack[k]);
2052 continue;
2053 }
2054 }
2055 // See whether there are parallel sides, e.g. a trapezium (UK English).
2056 if (np >= 4 && SplitTrapezium(stack[k], stack, panelsOut, epsang)) continue;
2057
2058 // Find a right-angled corner we can cut off.
2059 if (m_debug) std::cout << " Trying to find a right-angle\n";
2060 bool corner = false;
2061 for (unsigned int ip = 0; ip < np; ++ip) {
2062 // Take only right angles.
2063 const unsigned int inext = NextPoint(ip, np);
2064 const unsigned int iprev = PrevPoint(ip, np);
2065
2066 const double dxprev = xp1[iprev] - xp1[ip];
2067 const double dyprev = yp1[iprev] - yp1[ip];
2068 const double dxnext = xp1[inext] - xp1[ip];
2069 const double dynext = yp1[inext] - yp1[ip];
2070 if (fabs(dxprev * dxnext + dyprev * dynext) >
2071 epsang * sqrt((dxprev * dxprev + dyprev * dyprev) *
2072 (dxnext * dxnext + dynext * dynext))) {
2073 continue;
2074 }
2075 // Ensure the midpoint is internal.
2076 if (np > 3) {
2077 const double xm = 0.5 * (xp1[iprev] + xp1[inext]);
2078 const double ym = 0.5 * (yp1[iprev] + yp1[inext]);
2079 bool inside = false, edge = false;
2080 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2081 if (!inside) continue;
2082 }
2083 // Check all vertex crossings.
2084 bool cross = false;
2085 for (unsigned int jp = 0; jp < np; ++jp) {
2086 // Accept immediate contact.
2087 const unsigned int jnext = NextPoint(jp, np);
2088 if (jp == iprev || jp == ip || jp == inext || jnext == iprev ||
2089 jnext == ip || jnext == inext)
2090 continue;
2091 // Check crossing.
2092 double xc = 0., yc = 0.;
2093 if (Crossing(xp1[iprev], yp1[iprev], xp1[inext], yp1[inext], xp1[jp],
2094 yp1[jp], xp1[jnext], yp1[jnext], xc, yc)) {
2095 cross = true;
2096 break;
2097 }
2098 }
2099 if (cross) continue;
2100 // Found a triangle, introduce shorthand node references.
2101 if (m_debug) {
2102 std::cout << " Cutting at right-angled corner " << ip << "\n";
2103 }
2104 corner = true;
2105 Panel panel = stack[k];
2106 panel.xv = {xp1[iprev], xp1[ip], xp1[inext]};
2107 panel.yv = {yp1[iprev], yp1[ip], yp1[inext]};
2108 stack.push_back(std::move(panel));
2109 // Eliminate this node from the polygon.
2110 stack[k].xv.erase(stack[k].xv.begin() + ip);
2111 stack[k].yv.erase(stack[k].yv.begin() + ip);
2112 stack.push_back(std::move(stack[k]));
2113 if (m_debug) {
2114 std::cout << " Going for another pass, NP = " << np << "\n";
2115 }
2116 break;
2117 }
2118 if (corner) continue;
2119
2120 // Find any corner we can cut off.
2121 if (m_debug) std::cout << " Trying to find a corner\n";
2122 corner = false;
2123 for (unsigned int ip = 0; ip < np; ++ip) { // 20
2124 const unsigned int iprev = PrevPoint(ip, np);
2125 const unsigned int inext = NextPoint(ip, np);
2126 // Ensure the midpoint is internal.
2127 if (np > 3) {
2128 const double xm = 0.5 * (xp1[iprev] + xp1[inext]);
2129 const double ym = 0.5 * (yp1[iprev] + yp1[inext]);
2130 bool inside = false, edge = false;
2131 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2132 if (!inside) continue;
2133 }
2134 // Check all vertex crossings.
2135 bool cross = false;
2136 for (unsigned int jp = 0; jp < np; ++jp) {
2137 const unsigned int jj = NextPoint(jp, np);
2138 // Accept immediate contact.
2139 if (jp == iprev || jp == ip || jp == inext || jj == iprev || jj == ip ||
2140 jj == inext)
2141 continue;
2142 // Check crossing.
2143 double xc = 0., yc = 0.;
2144 if (Crossing(xp1[iprev], yp1[iprev], xp1[inext], yp1[inext], xp1[jp],
2145 yp1[jp], xp1[jj], yp1[jj], xc, yc)) {
2146 cross = true;
2147 break;
2148 }
2149 }
2150 if (cross) continue;
2151 // Found a triangle, introduce shorthand node references.
2152 if (m_debug) std::cout << " Cutting at corner " << ip << "\n";
2153 corner = true;
2154 const double x1 = xp1[iprev];
2155 const double x2 = xp1[ip];
2156 const double x3 = xp1[inext];
2157 const double y1 = yp1[iprev];
2158 const double y2 = yp1[ip];
2159 const double y3 = yp1[inext];
2160 const double s21 = (x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1);
2161 const double s31 = (x3 - x1) * (x3 - x1) + (y3 - y1) * (y3 - y1);
2162 const double s32 = (x3 - x2) * (x3 - x2) + (y3 - y2) * (y3 - y2);
2163 const double s12 = (x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2);
2164 const double s13 = (x1 - x3) * (x1 - x3) + (y1 - y3) * (y1 - y3);
2165 const double s23 = (x2 - x3) * (x2 - x3) + (y2 - y3) * (y2 - y3);
2166 // Find the biggest opening angle.
2167 const double a1 =
2168 ((x2 - x1) * (x3 - x1) + (y2 - y1) * (y3 - y1)) / sqrt(s21 * s31);
2169 const double a2 =
2170 ((x3 - x2) * (x1 - x2) + (y3 - y2) * (y1 - y2)) / sqrt(s32 * s12);
2171 const double a3 =
2172 ((x1 - x3) * (x2 - x3) + (y1 - y3) * (y2 - y3)) / sqrt(s13 * s23);
2173 if (m_debug) {
2174 const double phi1 = acos(a1);
2175 const double phi2 = acos(a2);
2176 const double phi3 = acos(a3);
2177 std::cout << " Angles: " << RadToDegree * phi1 << ", "
2178 << RadToDegree * phi2 << ", " << RadToDegree * phi3 << "\n";
2179 const double sumphi = phi1 + phi2 + phi3;
2180 std::cout << " Sum = " << RadToDegree * sumphi << "\n";
2181 }
2182 // See whether one angle is more or less right-angled.
2183 if (fabs(a1) < epsang || fabs(a2) < epsang || fabs(a3) < epsang) {
2184 if (m_debug) std::cout << " Right-angled corner cut off.\n";
2185 Panel panel = stack[k];
2186 if (fabs(a1) < epsang) {
2187 panel.xv = {x3, x1, x2};
2188 panel.yv = {y3, y1, y2};
2189 } else if (fabs(a2) < epsang) {
2190 panel.xv = {x1, x2, x3};
2191 panel.yv = {y1, y2, y3};
2192 } else {
2193 panel.xv = {x2, x3, x1};
2194 panel.yv = {y2, y3, y1};
2195 }
2196 stack.push_back(std::move(panel));
2197 } else if (a1 <= a2 && a1 <= a3) {
2198 if (m_debug) std::cout << " A1 < A2, A3 - adding 2 triangles.\n";
2199 const double xc = x2 + a2 * (x3 - x2) * sqrt(s12 / s32);
2200 const double yc = y2 + a2 * (y3 - y2) * sqrt(s12 / s32);
2201 Panel panel1 = stack[k];
2202 panel1.xv = {x3, xc, x1};
2203 panel1.yv = {y3, yc, y1};
2204 stack.push_back(std::move(panel1));
2205 Panel panel2 = stack[k];
2206 panel2.xv = {x2, xc, x1};
2207 panel2.yv = {y2, yc, y1};
2208 stack.push_back(std::move(panel2));
2209 } else if (a2 <= a1 && a2 <= a3) {
2210 if (m_debug) std::cout << " A2 < A1, A3 - adding 2 triangles.\n";
2211 const double xc = x3 + a3 * (x1 - x3) * sqrt(s23 / s13);
2212 const double yc = y3 + a3 * (y1 - y3) * sqrt(s23 / s13);
2213 Panel panel1 = stack[k];
2214 panel1.xv = {x1, xc, x2};
2215 panel1.yv = {y1, yc, y2};
2216 stack.push_back(std::move(panel1));
2217 Panel panel2 = stack[k];
2218 panel2.xv = {x3, xc, x2};
2219 panel2.yv = {y3, yc, y2};
2220 stack.push_back(std::move(panel2));
2221 } else {
2222 if (m_debug) std::cout << " A3 < A1, A2 - adding 2 triangles.\n";
2223 const double xc = x1 + a1 * (x2 - x1) * sqrt(s31 / s21);
2224 const double yc = y1 + a1 * (y2 - y1) * sqrt(s31 / s21);
2225 Panel panel1 = stack[k];
2226 panel1.xv = {x1, xc, x3};
2227 panel1.yv = {y1, yc, y3};
2228 stack.push_back(std::move(panel1));
2229 Panel panel2 = stack[k];
2230 panel2.xv = {x2, xc, x3};
2231 panel2.yv = {y2, yc, y3};
2232 stack.push_back(std::move(panel2));
2233 }
2234 // Eliminate this node from the polygon.
2235 stack[k].xv.erase(stack[k].xv.begin() + ip);
2236 stack[k].yv.erase(stack[k].yv.begin() + ip);
2237 stack.push_back(std::move(stack[k]));
2238 if (m_debug) {
2239 std::cout << " Going for another pass, NP = " << np << "\n";
2240 }
2241 break;
2242 }
2243 if (corner) continue;
2244 std::cerr << m_className << "::MakePrimitives:\n "
2245 << "Unable to identify a corner to cut, probably a degenerate "
2246 "polygon.\n";
2247 // Next stack element.
2248 }
2249 for (auto& panel : panelsOut) {
2250 panel.zv.assign(panel.xv.size(), zmean);
2251 }
2252 return true;
2253}
2254
2255bool ComponentNeBem3d::SplitTrapezium(const Panel panelIn,
2256 std::vector<Panel>& stack,
2257 std::vector<Panel>& panelsOut,
2258 const double epsang) const {
2259 const auto xp1 = panelIn.xv;
2260 const auto yp1 = panelIn.yv;
2261 const unsigned int np = xp1.size();
2262 for (unsigned int ip = 0; ip < np; ++ip) {
2263 const unsigned int inext = NextPoint(ip, np);
2264 const double xi0 = xp1[ip];
2265 const double yi0 = yp1[ip];
2266 const double xi1 = xp1[inext];
2267 const double yi1 = yp1[inext];
2268 const double dxi = xi0 - xi1;
2269 const double dyi = yi0 - yi1;
2270 const double si2 = dxi * dxi + dyi * dyi;
2271 for (unsigned int jp = ip + 2; jp < np; ++jp) {
2272 const unsigned int jnext = NextPoint(jp, np);
2273 // Skip adjacent segments.
2274 if (ip == jp || ip == jnext || inext == jp || inext == jnext) {
2275 continue;
2276 }
2277 const double xj0 = xp1[jp];
2278 const double yj0 = yp1[jp];
2279 const double xj1 = xp1[jnext];
2280 const double yj1 = yp1[jnext];
2281 const double dxj = xj0 - xj1;
2282 const double dyj = yj0 - yj1;
2283 const double sj2 = dxj * dxj + dyj * dyj;
2284 // Require parallelism.
2285 const double sij = sqrt(si2 * sj2);
2286 if (fabs(dxi * dxj + dyi * dyj + sij) > epsang * sij) continue;
2287 if (m_debug) {
2288 std::cout << " Found parallel sections: "
2289 << ip << ", " << jp << "\n";
2290 }
2291 // Avoid division by zero.
2292 if (sj2 <= 0 || si2 <= 0) {
2293 std::cerr << m_className << "::SplitTrapezium:\n "
2294 << "Zero norm segment found; skipped.\n";
2295 continue;
2296 }
2297 // Establish the cutting lines.
2298 const double xl1 =
2299 ((xi0 - xj0) * (xj1 - xj0) + (yi0 - yj0) * (yj1 - yj0)) / sj2;
2300 const double xl2 =
2301 ((xi1 - xj0) * (xj1 - xj0) + (yi1 - yj0) * (yj1 - yj0)) / sj2;
2302 const double xl3 =
2303 ((xj0 - xi0) * (xi1 - xi0) + (yj0 - yi0) * (yi1 - yi0)) / si2;
2304 const double xl4 =
2305 ((xj1 - xi0) * (xi1 - xi0) + (yj1 - yi0) * (yi1 - yi0)) / si2;
2306 if (m_debug) {
2307 std::cout << " xl1 = " << xl1 << ", xl2 = " << xl2 << ", "
2308 << "xl3 = " << xl3 << ", xl4 = " << xl4 << "\n";
2309 }
2310 // Check that there is at all a rectangle.
2311 const double r1 = (xl1 + epsang) * (1. + epsang - xl1);
2312 const double r2 = (xl2 + epsang) * (1. + epsang - xl2);
2313 const double r3 = (xl3 + epsang) * (1. + epsang - xl3);
2314 const double r4 = (xl4 + epsang) * (1. + epsang - xl4);
2315 if ((r1 < 0 && r4 < 0) || (r2 < 0 && r3 < 0)) {
2316 if (m_debug) std::cout << " No rectangle.\n";
2317 continue;
2318 }
2319 // Determine the rectangular part.
2320 std::vector<double> xpl(4, 0.);
2321 std::vector<double> ypl(4, 0.);
2322 if (r1 >= 0) {
2323 xpl[0] = xi0;
2324 ypl[0] = yi0;
2325 xpl[1] = xj0 + xl1 * (xj1 - xj0);
2326 ypl[1] = yj0 + xl1 * (yj1 - yj0);
2327 } else if (r4 >= 0) {
2328 xpl[0] = xi0 + xl4 * (xi1 - xi0);
2329 ypl[0] = yi0 + xl4 * (yi1 - yi0);
2330 xpl[1] = xj1;
2331 ypl[1] = yj1;
2332 }
2333 if (r2 >= 0) {
2334 xpl[2] = xj0 + xl2 * (xj1 - xj0);
2335 ypl[2] = yj0 + xl2 * (yj1 - yj0);
2336 xpl[3] = xi1;
2337 ypl[3] = yi1;
2338 } else if (r3 >= 0) {
2339 xpl[2] = xj0;
2340 ypl[2] = yj0;
2341 xpl[3] = xi0 + xl3 * (xi1 - xi0);
2342 ypl[3] = yi0 + xl3 * (yi1 - yi0);
2343 }
2344 // Verify that the midpoints of these lines are internal.
2345 double xm = 0.5 * (xpl[0] + xpl[1]);
2346 double ym = 0.5 * (ypl[0] + ypl[1]);
2347 bool inside = false, edge = false;
2348 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2349 if (!(inside || edge)) {
2350 if (m_debug) std::cout << " Midpoint 1 not internal.\n";
2351 continue;
2352 }
2353 xm = 0.5 * (xpl[2] + xpl[3]);
2354 ym = 0.5 * (ypl[2] + ypl[3]);
2355 Polygon::Inside(xp1, yp1, xm, ym, inside, edge);
2356 if (!(inside || edge)) {
2357 if (m_debug) std::cout << " Midpoint 2 not internal.\n";
2358 continue;
2359 }
2360
2361 const unsigned int iprev = PrevPoint(ip, np);
2362 const unsigned int jprev = PrevPoint(jp, np);
2363 // Ensure there are no crossings, accepting contact.
2364 bool cross = false;
2365 for (unsigned int i = 0; i < np; ++i) {
2366 if ((i == iprev && r1 >= 0) || i == ip || (i == inext && r2 >= 0) ||
2367 (i == jprev && r3 >= 0) || i == jp || (i == jnext && r4 >= 0))
2368 continue;
2369 const unsigned int ii = NextPoint(i, np);
2370 double xc = 0., yc = 0.;
2371 if (Crossing(xp1[i], yp1[i], xp1[ii], yp1[ii], xpl[0], ypl[0], xpl[1],
2372 ypl[1], xc, yc) ||
2373 Crossing(xp1[i], yp1[i], xp1[ii], yp1[ii], xpl[2], ypl[2], xpl[3],
2374 ypl[3], xc, yc)) {
2375 if (m_debug) {
2376 std::cout << " Crossing (edge " << i << ", " << ii
2377 << "), ip = " << ip << ", jp = " << jp << "\n";
2378 }
2379 cross = true;
2380 break;
2381 }
2382 }
2383 if (cross) continue;
2384 // Add the rectangular part.
2385 if ((fabs(xl1) < epsang && fabs(xl3) < epsang) ||
2386 (fabs(1. - xl2) < epsang && fabs(1. - xl4) < epsang)) {
2387 if (m_debug) std::cout << " Not stored, degenerate.\n";
2388 } else {
2389 if (m_debug) std::cout << " Adding rectangle.\n";
2390 Panel panel = panelIn;
2391 panel.xv = xpl;
2392 panel.yv = ypl;
2393 panelsOut.push_back(std::move(panel));
2394 }
2395 // First non-rectangular section.
2396 xpl.clear();
2397 ypl.clear();
2398 if (m_debug) std::cout << " First non-rectangular section.\n";
2399 for (unsigned int i = jp + 1; i <= ip + np; ++i) {
2400 const unsigned int ii = i % np;
2401 xpl.push_back(xp1[ii]);
2402 ypl.push_back(yp1[ii]);
2403 }
2404 if (r1 >= 0 && r4 >= 0) {
2405 if (m_debug) std::cout << " 1-4 degenerate\n";
2406 } else if (r1 >= 0) {
2407 if (m_debug) std::cout << " Using 1\n";
2408 xpl.push_back(xj0 + xl1 * (xj1 - xj0));
2409 ypl.push_back(yj0 + xl1 * (yj1 - yj0));
2410 } else if (r4 >= 0) {
2411 if (m_debug) std::cout << " Using 4\n";
2412 xpl.push_back(xi0 + xl4 * (xi1 - xi0));
2413 ypl.push_back(yi0 + xl4 * (yi1 - yi0));
2414 } else {
2415 if (m_debug) std::cout << " Neither 1 nor 4, should not happen\n";
2416 }
2417 if (xpl.size() < 3) {
2418 if (m_debug) {
2419 std::cout << " Not stored, only " << xpl.size()
2420 << " vertices.\n";
2421 }
2422 } else {
2423 if (m_debug) {
2424 std::cout << " Adding non-rectangular part with " << xpl.size()
2425 << " vertices.\n";
2426 }
2427 Panel panel = panelIn;
2428 panel.xv = xpl;
2429 panel.yv = ypl;
2430 stack.push_back(std::move(panel));
2431 }
2432 // Second non-rectangular section.
2433 xpl.clear();
2434 ypl.clear();
2435 if (m_debug) std::cout << " Second non-rectangular section.\n";
2436 for (unsigned int i = ip + 1; i <= jp; ++i) {
2437 const unsigned int ii = i % np;
2438 xpl.push_back(xp1[ii]);
2439 ypl.push_back(yp1[ii]);
2440 }
2441 if (r2 >= 0 && r3 >= 0) {
2442 if (m_debug) std::cout << " 2-3 degenerate\n";
2443 } else if (r2 >= 0) {
2444 if (m_debug) std::cout << " Using 2\n";
2445 xpl.push_back(xj0 + xl2 * (xj1 - xj0));
2446 ypl.push_back(yj0 + xl2 * (yj1 - yj0));
2447 } else if (r3 >= 0) {
2448 if (m_debug) std::cout << " Using 3\n";
2449 xpl.push_back(xi0 + xl3 * (xi1 - xi0));
2450 ypl.push_back(yi0 + xl3 * (yi1 - yi0));
2451 } else {
2452 if (m_debug) std::cout << " Neither 2 nor 3, should not happen\n";
2453 }
2454 if (xpl.size() < 3) {
2455 if (m_debug) {
2456 std::cout << " Not stored, only " << xpl.size()
2457 << " vertices.\n";
2458 }
2459 } else {
2460 if (m_debug) {
2461 std::cout << " Adding non-rectangular part with " << xpl.size()
2462 << " vertices.\n";
2463 }
2464 Panel panel = panelIn;
2465 panel.xv = xpl;
2466 panel.yv = ypl;
2467 stack.push_back(std::move(panel));
2468 }
2469 return true;
2470 }
2471 }
2472 return false;
2473}
2474
2475bool ComponentNeBem3d::GetPrimitive(const unsigned int i, double& a, double& b,
2476 double& c, std::vector<double>& xv,
2477 std::vector<double>& yv,
2478 std::vector<double>& zv, int& interface,
2479 double& v, double& q,
2480 double& lambda) const {
2481 if (i >= m_primitives.size()) {
2482 std::cerr << m_className << "::GetPrimitive: Index out of range.\n";
2483 return false;
2484 }
2485 const auto& primitive = m_primitives[i];
2486 a = primitive.a;
2487 b = primitive.b;
2488 c = primitive.c;
2489 xv = primitive.xv;
2490 yv = primitive.yv;
2491 zv = primitive.zv;
2492 interface = primitive.interface;
2493 v = primitive.v;
2494 q = primitive.q;
2495 lambda = primitive.lambda;
2496 return true;
2497}
2498
2499bool ComponentNeBem3d::GetElement(const unsigned int i,
2500 std::vector<double>& xv,
2501 std::vector<double>& yv,
2502 std::vector<double>& zv, int& interface,
2503 double& bc, double& lambda) const {
2504
2505 if (i >= m_elements.size()) {
2506 std::cerr << m_className << "::GetElement: Index out of range.\n";
2507 return false;
2508 }
2509 const auto& element = m_elements[i];
2510 xv = element.xv;
2511 yv = element.yv;
2512 zv = element.zv;
2513 interface = element.interface;
2514 bc = element.bc;
2515 lambda = element.lambda;
2516 return true;
2517}
2518
2520 m_primitives.clear();
2521 m_elements.clear();
2522 m_ready = false;
2523}
2524
2526 std::cerr << m_className << "::UpdatePeriodicity:\n"
2527 << " Periodicities are not supported.\n";
2528}
2529}
Abstract base class for components.
GeometryBase * m_geometry
Pointer to the geometry.
bool m_ready
Ready for use?
std::string m_className
Class name.
bool m_debug
Switch on/off debugging messages.
virtual Medium * GetMedium(const double x, const double y, const double z)
Get the medium at a given location (x, y, z).
void Reset() override
Reset the component.
void UpdatePeriodicity() override
Verify periodicities.
bool GetVoltageRange(double &vmin, double &vmax) override
Calculate the voltage range [V].
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)
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 ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
virtual unsigned int GetNumberOfSolids() const
Return the number of solids in the geometry.
Definition: GeometryBase.hh:25
virtual Solid * GetSolid(const unsigned int) const
Get a solid from the list.
Definition: GeometryBase.hh:27
Abstract base class for media.
Definition: Medium.hh:13
double GetDielectricConstant() const
Get the relative static dielectric constant.
Definition: Medium.hh:42
@ DielectricCharge
Definition: Solid.hh:130
@ PerpendicularField
Definition: Solid.hh:132
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:135
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
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
std::vector< double > xv
X-coordinates of vertices.
Definition: Solid.hh:15
std::vector< double > yv
Y-coordinates of vertices.
Definition: Solid.hh:17