OpenLB 1.8.1
Loading...
Searching...
No Matches
STLreaderForSubgridParticleWallContact.hh
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2010-2024 Thomas Henn, Mathias J. Krause, Jonathan Jeppener-Haltenhoff, Christoph Gaul
4 * E-mail contact: info@openlb.net
5 * The most recent release of OpenLB can be downloaded at
6 * <http://www.openlb.net/>
7 *
8 * This program is free software; you can redistribute it and/or
9 * modify it under the terms of the GNU General Public License
10 * as published by the Free Software Foundation; either version 2
11 * of the License, or (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public
19 * License along with this program; if not, write to the Free
20 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
21 * Boston, MA 02110-1301, USA.
22 */
23
28#ifndef STL_READER_HH
29#define STL_READER_HH
30
31
32#include <iostream>
33#include <sstream>
34#include <fstream>
35#include <stdexcept>
36#include "core/singleton.h"
38#include "octree.hh"
39#include "stlReader.h"
40
41#ifdef FEATURE_VTK
42#include <vtkSmartPointer.h>
43#include <vtkUnstructuredGrid.h>
44#include <vtkPoints.h>
45#include <vtkCellArray.h>
46#include <vtkFloatArray.h>
47#include <vtkPointData.h>
48#include <vtkTriangle.h>
49#include <vtkXMLUnstructuredGridWriter.h>
50#endif
51
52// All OpenLB code is contained in this namespace.
53namespace olb {
54
55template<typename T>
57{
58 Vector<T,3> A = point[0];
59 Vector<T,3> B = point[1];
60 Vector<T,3> C = point[2];
61 Vector<T,3> b, c;
62 T bb = 0., bc = 0., cc = 0.;
63
64 for (int i = 0; i < 3; i++) {
65 b[i] = B[i] - A[i];
66 c[i] = C[i] - A[i];
67 bb += b[i] * b[i];
68 bc += b[i] * c[i];
69 cc += c[i] * c[i];
70 }
71
72 normal[0] = b[1] * c[2] - b[2] * c[1];
73 normal[1] = b[2] * c[0] - b[0] * c[2];
74 normal[2] = b[0] * c[1] - b[1] * c[0];
75
76 T norm = util::sqrt(
77 util::pow(normal[0], 2) + util::pow(normal[1], 2) + util::pow(normal[2], 2));
78 normal[0] /= norm;
79 normal[1] /= norm;
80 normal[2] /= norm;
81
82 T D = 1.0 / (cc * bb - bc * bc);
83 T bbD = bb * D;
84 T bcD = bc * D;
85 T ccD = cc * D;
86
87 kBeta = 0.;
88 kGamma = 0.;
89 d = 0.;
90
91 for (int i = 0; i < 3; i++) {
92 uBeta[i] = b[i] * ccD - c[i] * bcD;
93 uGamma[i] = c[i] * bbD - b[i] * bcD;
94 kBeta -= A[i] * uBeta[i];
95 kGamma -= A[i] * uGamma[i];
96 d += A[i] * normal[i];
97 }
98}
99
100template<typename T>
101Vector<T,3> STLtriangle<T>::getCenter()
102{
103 Vector<T,3> center( T(0) );
104
105 center[0] = (point[0][0] + point[1][0]
106 + point[2][0]) / 3.;
107 center[1] = (point[0][1] + point[1][1]
108 + point[2][1]) / 3.;
109 center[2] = (point[0][2] + point[1][2]
110 + point[2][2]) / 3.;
111
112 return center;
113}
114
115template<typename T>
116std::vector<T> STLtriangle<T>::getE0()
117{
118 Vector<T,3> vec;
119 vec[0] = point[0][0] - point[1][0];
120 vec[1] = point[0][1] - point[1][1];
121 vec[2] = point[0][2] - point[1][2];
122 return vec;
123}
124
125template<typename T>
126std::vector<T> STLtriangle<T>::getE1()
127{
128 Vector<T,3> vec;
129 vec[0] = point[0][0] - point[2][0];
130 vec[1] = point[0][1] - point[2][1];
131 vec[2] = point[0][2] - point[2][2];
132 return vec;
133}
134
135template<typename T>
136bool STLtriangle<T>::isPointInside(const PhysR<T,3>& pt) const
137{
138 // tests with T=double and T=float show that the epsilon must be increased
139 const T epsilon = std::numeric_limits<BaseType<T>>::epsilon()*T(10);
140
141 const T beta = pt * uBeta + kBeta;
142 const T gamma = pt * uGamma + kGamma;
143
144 // check if approximately equal
145 if ( util::nearZero(norm(pt - (point[0] + beta*(point[1]-point[0]) + gamma*(point[2]-point[0]))), epsilon) ) {
146 const T alpha = T(1) - beta - gamma;
147 return (beta >= T(0) || util::nearZero(beta, epsilon))
148 && (gamma >= T(0) || util::nearZero(gamma, epsilon))
149 && (alpha >= T(0) || util::nearZero(alpha, epsilon));
150 }
151 return false;
152}
153
154/* Schnitttest nach
155 * http://www.uninformativ.de/bin/RaytracingSchnitttests-76a577a-CC-BY.pdf
156 *
157 * Creative Commons Namensnennung 3.0 Deutschland
158 * http://creativecommons.org/licenses/by/3.0/de/
159 *
160 * P. Hofmann, 22. August 2010
161 *
162 */
163template<typename T>
164bool STLtriangle<T>::testRayIntersect(const Vector<T,3>& pt,
165 const Vector<T,3>& dir,
166 Vector<T,3>& q, T& alpha, const T& rad,
167 bool print)
168{
169 T rn = 0.;
170 Vector<T,3> testPt = pt + rad * normal;
171 /* Vector<T,3> help; */
172
173 for (int i = 0; i < 3; i++) {
174 rn += dir[i] * normal[i];
175 }
176#ifdef OLB_DEBUG
177
178 if (print) {
179 std::cout << "Pt: " << pt[0] << " " << pt[1] << " " << pt[2] << std::endl;
180 }
181 if (print)
182 std::cout << "testPt: " << testPt[0] << " " << testPt[1] << " " << testPt[2]
183 << std::endl;
184 if (print)
185 std::cout << "PosNeg: "
186 << normal[0] * testPt[0] + normal[1] * testPt[1] + normal[2] * testPt[2]
187 - d << std::endl;
188 if (print)
189 std::cout << "Normal: " << normal[0] << " " << normal[1] << " " << normal[2]
190 << std::endl;
191#endif
192
193 // Schnitttest Flugrichtung -> Ebene
194 if (util::fabs(rn) < std::numeric_limits<T>::epsilon()) {
195#ifdef OLB_DEBUG
196 if (print) {
197 std::cout << "FALSE 1" << std::endl;
198 }
199#endif
200 return false;
201 }
202 alpha = d - testPt[0] * normal[0] - testPt[1] * normal[1] - testPt[2] * normal[2];
203 // alpha -= testPt[i] * normal[i];
204 alpha /= rn;
205
206 // Abstand Partikel Ebene
207 if (alpha < -std::numeric_limits<T>::epsilon()) {
208#ifdef OLB_DEBUG
209 if (print) {
210 std::cout << "FALSE 2" << std::endl;
211 }
212#endif
213 return false;
214 }
215 for (int i = 0; i < 3; i++) {
216 q[i] = testPt[i] + alpha * dir[i];
217 }
218 T beta = kBeta;
219 for (int i = 0; i < 3; i++) {
220 beta += uBeta[i] * q[i];
221 }
222#ifdef OLB_DEBUG
223 T dist = util::sqrt(
224 util::pow(q[0] - testPt[0], 2) + util::pow(q[1] - testPt[1], 2)
225 + util::pow(q[2] - testPt[2], 2));
226#endif
227
228 // Schnittpunkt q in der Ebene?
229 if (beta < -std::numeric_limits<T>::epsilon()) {
230#ifdef OLB_DEBUG
231
232 if (print) {
233 std::cout << "FALSE 3 BETA " << beta << " DIST " << dist << std::endl;
234 }
235#endif
236 return false;
237 }
238 T gamma = kGamma;
239 for (int i = 0; i < 3; i++) {
240 gamma += uGamma[i] * q[i];
241 }
242 if (gamma < -std::numeric_limits<T>::epsilon()) {
243#ifdef OLB_DEBUG
244 if (print) {
245 std::cout << "FALSE 4 GAMMA " << gamma << " DIST " << dist << std::endl;
246 }
247#endif
248 return false;
249 }
250 if (1. - beta - gamma < -std::numeric_limits<T>::epsilon()) {
251#ifdef OLB_DEBUG
252 if (print)
253 std::cout << "FALSE 5 VAL " << 1 - beta - gamma << " DIST " << dist
254 << std::endl;
255#endif
256 return false;
257 }
258#ifdef OLB_DEBUG
259 if (print) {
260 std::cout << "TRUE" << " GAMMA " << gamma << " BETA " << beta << std::endl;
261 }
262#endif
263 return true;
264}
265
270template<typename T>
272 const Vector<T,3>& pt) const
273{
274
275 const T nEps = -std::numeric_limits<T>::epsilon();
276 const T Eps = std::numeric_limits<T>::epsilon();
277
278 Vector<T,3> ab = point[1] - point[0];
279 Vector<T,3> ac = point[2] - point[0];
280 Vector<T,3> bc = point[2] - point[1];
281
282 T snom = (pt - point[0])*ab;
283 T sdenom = (pt - point[1])*(point[0] - point[1]);
284
285 T tnom = (pt - point[0])*ac;
286 T tdenom = (pt - point[2])*(point[0] - point[2]);
287
288 if (snom < nEps && tnom < nEps) {
289 return point[0];
290 }
291
292 T unom = (pt - point[1])*bc;
293 T udenom = (pt - point[2])*(point[1] - point[2]);
294
295 if (sdenom < nEps && unom < nEps) {
296 return point[1];
297 }
298 if (tdenom < nEps && udenom < nEps) {
299 return point[2];
300 }
301
302 T vc = normal*crossProduct3D(point[0] - pt, point[1] - pt);
303
304 if (vc < nEps && snom > Eps && sdenom > Eps) {
305 return point[0] + snom / (snom + sdenom) * ab;
306 }
307
308 T va = normal*crossProduct3D(point[1] - pt, point[2] - pt);
309
310 if (va < nEps && unom > Eps && udenom > Eps) {
311 return point[1] + unom / (unom + udenom) * bc;
312 }
313
314 T vb = normal*crossProduct3D(point[2] - pt, point[0] - pt);
315
316 if (vb < nEps && tnom > Eps && tdenom > Eps) {
317 return point[0] + tnom / (tnom + tdenom) * ac;
318 }
319
320 T u = va / (va + vb + vc);
321 T v = vb / (va + vb + vc);
322 T w = 1. - u - v;
323
324 return u * point[0] + v * point[1] + w * point[2];
325}
326
327template<typename T>
328bool STLtriangle<T>::getPointToEdgeDistances(const Vector<T,3>& input, Vector<T,3>& output, T sensitivity)
329{
330 auto P1 = point[0];
331 auto P2 = point[1];
332 output[0] = norm(crossProduct3D(P2-input,P1-P2))/norm(P1-P2);
333 P1 = point[0];
334 P2 = point[2];
335 output[1] = norm(crossProduct3D(P2-input,P1-P2))/norm(P1-P2);
336 P1 = point[1];
337 P2 = point[2];
338 output[2] = norm(crossProduct3D(P2-input,P1-P2))/norm(P1-P2);
339 return true;
340}
341
342 template<typename T>
343 bool STLtriangle<T>::isEdgePoint (const Vector<T,3> & input, Vector<T,3>& P1, Vector<T,3> & P2, T sensitivity )
344 {
345 using namespace std;
346 if(input[0] < sensitivity && input[1] >= sensitivity && input[2] >= sensitivity )
347 {
348 P1 = point[0];
349 P2 = point[1];
350 //cout << "edge point " << P1 << " " << P2 << std::endl;
351 return true;
352 }
353 if(input[0] >= sensitivity && input[1] < sensitivity && input[2] >= sensitivity )
354 {
355 P1 = point[0];
356 P2 = point[2];
357 //cout << "edge point" << P1 << " " << P2 << std::endl;
358 return true;
359 }
360 if(input[0] >= sensitivity && input[1] >= sensitivity && input[2] < sensitivity )
361 {
362 P1 = point[1];
363 P2 = point[2];
364 //cout << "edge point" << P1 << " " << P2 << std::endl;
365 return true;
366 }
367 return false;
368 }
369
370 template<typename T>
371 bool STLtriangle<T>::isVortexPoint (const Vector<T,3> & input,Vector<T,3>& P, T sensitivity )
372 {
373 olb::OstreamManager clout = OstreamManager(std::cout, "STLtriangle");
374 using namespace std;
375 if (input[0] < sensitivity && input[1] < sensitivity && input[2] < sensitivity )
376 {
377 cout << "Error isVortexPoint! Possible reduce sensitivity!" << std::endl;
378 return false;
379 }
380 if(input[0] < sensitivity && input[1] < sensitivity)
381 {
382 P = point[0];
383 //cout << "vortex point " << P << std::endl;
384 return true;
385 }
386 if(input[0] < sensitivity && input[2] < sensitivity)
387 {
388 P = point[1];
389 //cout << "vortex point " << P << std::endl;
390 return true;
391 }
392 if(input[1] < sensitivity && input[2] < sensitivity)
393 {
394 P = point[2];
395 //cout << "vortex point " << P << std::endl;
396 return true;
397 }
398 return false;
399 }
400
401
402template<typename T>
403STLmesh<T>::STLmesh(std::string fName, T stlSize)
404 : _fName(fName),
405 _min(T()),
406 _max(T()),
407 _maxDist2(0),
408 clout(std::cout, "STLmesh")
409{
410 std::ifstream f(fName.c_str(), std::ios::in);
411 _triangles.reserve(10000);
412 if (!f.good()) {
413 throw std::runtime_error("STL File not valid.");
414 }
415 char buf[6];
416 buf[5] = 0;
417 f.read(buf, 5);
418 const std::string asciiHeader = "solid";
419 if (std::string(buf) == asciiHeader) {
420 f.seekg(0, std::ios::beg);
421 if (f.good()) {
422 std::string s0, s1;
423 int i = 0;
424 while (!f.eof()) {
425 f >> s0;
426 if (s0 == "facet") {
427 STLtriangle<T> tri;
428 f >> s1 >> tri.normal[0] >> tri.normal[1] >> tri.normal[2];
429 f >> s0 >> s1;
430 f >> s0 >> tri.point[0][0] >> tri.point[0][1]
431 >> tri.point[0][2];
432 f >> s0 >> tri.point[1][0] >> tri.point[1][1]
433 >> tri.point[1][2];
434 f >> s0 >> tri.point[2][0] >> tri.point[2][1]
435 >> tri.point[2][2];
436 f >> s0;
437 f >> s0;
438 for (int k = 0; k < 3; k++) {
439 tri.point[0][k] *= stlSize;
440 tri.point[1][k] *= stlSize;
441 tri.point[2][k] *= stlSize;
442 }
443 if (i == 0) {
444 _min = T();
445 _max = T();
446
447 _min[0] = tri.point[0][0];
448 _min[1] = tri.point[0][1];
449 _min[2] = tri.point[0][2];
450
451 _max[0] = tri.point[0][0];
452 _max[1] = tri.point[0][1];
453 _max[2] = tri.point[0][2];
454
455 _min[0] = util::min(_min[0], tri.point[1][0]);
456 _min[1] = util::min(_min[1], tri.point[1][1]);
457 _min[2] = util::min(_min[2], tri.point[1][2]);
458
459 _max[0] = util::max(_max[0], tri.point[1][0]);
460 _max[1] = util::max(_max[1], tri.point[1][1]);
461 _max[2] = util::max(_max[2], tri.point[1][2]);
462
463 _min[0] = util::min(_min[0], tri.point[2][0]);
464 _min[1] = util::min(_min[1], tri.point[2][1]);
465 _min[2] = util::min(_min[2], tri.point[2][2]);
466
467 _max[0] = util::max(_max[0], tri.point[2][0]);
468 _max[1] = util::max(_max[1], tri.point[2][1]);
469 _max[2] = util::max(_max[2], tri.point[2][2]);
470
471 }
472 else {
473 _min[0] = util::min(_min[0], tri.point[0][0]);
474 _min[1] = util::min(_min[1], tri.point[0][1]);
475 _min[2] = util::min(_min[2], tri.point[0][2]);
476
477 _max[0] = util::max(_max[0], tri.point[0][0]);
478 _max[1] = util::max(_max[1], tri.point[0][1]);
479 _max[2] = util::max(_max[2], tri.point[0][2]);
480
481 _min[0] = util::min(_min[0], tri.point[1][0]);
482 _min[1] = util::min(_min[1], tri.point[1][1]);
483 _min[2] = util::min(_min[2], tri.point[1][2]);
484
485 _max[0] = util::max(_max[0], tri.point[1][0]);
486 _max[1] = util::max(_max[1], tri.point[1][1]);
487 _max[2] = util::max(_max[2], tri.point[1][2]);
488
489 _min[0] = util::min(_min[0], tri.point[2][0]);
490 _min[1] = util::min(_min[1], tri.point[2][1]);
491 _min[2] = util::min(_min[2], tri.point[2][2]);
492
493 _max[0] = util::max(_max[0], tri.point[2][0]);
494 _max[1] = util::max(_max[1], tri.point[2][1]);
495 _max[2] = util::max(_max[2], tri.point[2][2]);
496 }
497
498 i++;
499 tri.init();
500 _triangles.push_back(tri);
501
502 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[1]),
503 _maxDist2);
504 _maxDist2 = util::max(distPoints(tri.point[2], tri.point[1]),
505 _maxDist2);
506 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[2]),
507 _maxDist2);
508 }
509 else if (s0 == "endsolid") {
510 break;
511 }
512 }
513 }
514 }
515 else {
516 f.close();
517 f.open(fName.c_str(), std::ios::in | std::ios::binary);
518 char comment[80];
519 f.read(comment, 80);
520
521 if (!f.good()) {
522 throw std::runtime_error("STL File not valid.");
523 }
524
525 comment[79] = 0;
526 int32_t nFacets;
527 f.read(reinterpret_cast<char *>(&nFacets), sizeof(int32_t));
528
529 if (!f.good()) {
530 throw std::runtime_error("STL File not valid.");
531 }
532
533 float v[12];
534 std::uint16_t uint16;
535 for (int32_t i = 0; i < nFacets; ++i) {
536 for (unsigned int j = 0; j < 12; ++j) {
537 f.read(reinterpret_cast<char *>(&v[j]), sizeof(float));
538 }
539 f.read(reinterpret_cast<char *>(&uint16), sizeof(std::uint16_t));
540 STLtriangle<T> tri;
541 tri.normal[0] = v[0];
542 tri.normal[1] = v[1];
543 tri.normal[2] = v[2];
544 tri.point[0][0] = v[3];
545 tri.point[0][1] = v[4];
546 tri.point[0][2] = v[5];
547 tri.point[1][0] = v[6];
548 tri.point[1][1] = v[7];
549 tri.point[1][2] = v[8];
550 tri.point[2][0] = v[9];
551 tri.point[2][1] = v[10];
552 tri.point[2][2] = v[11];
553
554 for (int k = 0; k < 3; k++) {
555 tri.point[0][k] *= stlSize;
556 tri.point[1][k] *= stlSize;
557 tri.point[2][k] *= stlSize;
558 }
559 if (i == 0) {
560 _min[0] = tri.point[0][0];
561 _min[1] = tri.point[0][1];
562 _min[2] = tri.point[0][2];
563
564 _max[0] = tri.point[0][0];
565 _max[1] = tri.point[0][1];
566 _max[2] = tri.point[0][2];
567
568 _min[0] = util::min(_min[0], (T) tri.point[1][0]);
569 _min[1] = util::min(_min[1], (T) tri.point[1][1]);
570 _min[2] = util::min(_min[2], (T) tri.point[1][2]);
571
572 _max[0] = util::max(_max[0], (T) tri.point[1][0]);
573 _max[1] = util::max(_max[1], (T) tri.point[1][1]);
574 _max[2] = util::max(_max[2], (T) tri.point[1][2]);
575
576 _min[0] = util::min(_min[0], (T) tri.point[2][0]);
577 _min[1] = util::min(_min[1], (T) tri.point[2][1]);
578 _min[2] = util::min(_min[2], (T) tri.point[2][2]);
579
580 _max[0] = util::max(_max[0], (T) tri.point[2][0]);
581 _max[1] = util::max(_max[1], (T) tri.point[2][1]);
582 _max[2] = util::max(_max[2], (T) tri.point[2][2]);
583
584 }
585 else {
586 _min[0] = util::min(_min[0], (T) tri.point[0][0]);
587 _min[1] = util::min(_min[1], (T) tri.point[0][1]);
588 _min[2] = util::min(_min[2], (T) tri.point[0][2]);
589
590 _max[0] = util::max(_max[0], (T) tri.point[0][0]);
591 _max[1] = util::max(_max[1], (T) tri.point[0][1]);
592 _max[2] = util::max(_max[2], (T) tri.point[0][2]);
593
594 _min[0] = util::min(_min[0], (T) tri.point[1][0]);
595 _min[1] = util::min(_min[1], (T) tri.point[1][1]);
596 _min[2] = util::min(_min[2], (T) tri.point[1][2]);
597
598 _max[0] = util::max(_max[0], (T) tri.point[1][0]);
599 _max[1] = util::max(_max[1], (T) tri.point[1][1]);
600 _max[2] = util::max(_max[2], (T) tri.point[1][2]);
601
602 _min[0] = util::min(_min[0], (T) tri.point[2][0]);
603 _min[1] = util::min(_min[1], (T) tri.point[2][1]);
604 _min[2] = util::min(_min[2], (T) tri.point[2][2]);
605
606 _max[0] = util::max(_max[0], (T) tri.point[2][0]);
607 _max[1] = util::max(_max[1], (T) tri.point[2][1]);
608 _max[2] = util::max(_max[2], (T) tri.point[2][2]);
609 }
610 tri.init();
611 _triangles.push_back(tri);
612
613 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[1]), _maxDist2);
614 _maxDist2 = util::max(distPoints(tri.point[2], tri.point[1]), _maxDist2);
615 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[2]), _maxDist2);
616 }
617 }
618 f.close();
619}
620
621template<typename T>
622STLmesh<T>::STLmesh(const std::vector<std::vector<T>> meshPoints, T stlSize)
623 : _fName("meshPoints.stl"),
624 _min(T()),
625 _max(T()),
626 _maxDist2(0),
627 clout(std::cout, "STLmesh")
628{
629 _triangles.reserve(10000);
630 for (size_t i = 0; i < meshPoints.size() / 3; i++) {
631 STLtriangle<T> tri;
632 tri.point[0][0] = meshPoints[i*3 + 0][0];
633 tri.point[0][1] = meshPoints[i*3 + 0][1];
634 tri.point[0][2] = meshPoints[i*3 + 0][2];
635
636 tri.point[1][0] = meshPoints[i*3 + 1][0];
637 tri.point[1][1] = meshPoints[i*3 + 1][1];
638 tri.point[1][2] = meshPoints[i*3 + 1][2];
639
640 tri.point[2][0] = meshPoints[i*3 + 2][0];
641 tri.point[2][1] = meshPoints[i*3 + 2][1];
642 tri.point[2][2] = meshPoints[i*3 + 2][2];
643 for (int k = 0; k < 3; k++) {
644 tri.point[0][k] *= stlSize;
645 tri.point[1][k] *= stlSize;
646 tri.point[2][k] *= stlSize;
647 }
648 if (i == 0) {
649 _min*=T();
650 _max*=T();
651
652 _min[0] = tri.point[0][0];
653 _min[1] = tri.point[0][1];
654 _min[2] = tri.point[0][2];
655
656 _max[0] = tri.point[0][0];
657 _max[1] = tri.point[0][1];
658 _max[2] = tri.point[0][2];
659
660 _min[0] = util::min(_min[0], (T) tri.point[1][0]);
661 _min[1] = util::min(_min[1], (T) tri.point[1][1]);
662 _min[2] = util::min(_min[2], (T) tri.point[1][2]);
663
664 _max[0] = util::max(_max[0], (T) tri.point[1][0]);
665 _max[1] = util::max(_max[1], (T) tri.point[1][1]);
666 _max[2] = util::max(_max[2], (T) tri.point[1][2]);
667
668 _min[0] = util::min(_min[0], (T) tri.point[2][0]);
669 _min[1] = util::min(_min[1], (T) tri.point[2][1]);
670 _min[2] = util::min(_min[2], (T) tri.point[2][2]);
671
672 _max[0] = util::max(_max[0], (T) tri.point[2][0]);
673 _max[1] = util::max(_max[1], (T) tri.point[2][1]);
674 _max[2] = util::max(_max[2], (T) tri.point[2][2]);
675
676 }
677 else {
678 _min[0] = util::min(_min[0], (T) tri.point[0][0]);
679 _min[1] = util::min(_min[1], (T) tri.point[0][1]);
680 _min[2] = util::min(_min[2], (T) tri.point[0][2]);
681
682 _max[0] = util::max(_max[0], (T) tri.point[0][0]);
683 _max[1] = util::max(_max[1], (T) tri.point[0][1]);
684 _max[2] = util::max(_max[2], (T) tri.point[0][2]);
685
686 _min[0] = util::min(_min[0], (T) tri.point[1][0]);
687 _min[1] = util::min(_min[1], (T) tri.point[1][1]);
688 _min[2] = util::min(_min[2], (T) tri.point[1][2]);
689
690 _max[0] = util::max(_max[0], (T) tri.point[1][0]);
691 _max[1] = util::max(_max[1], (T) tri.point[1][1]);
692 _max[2] = util::max(_max[2], (T) tri.point[1][2]);
693
694 _min[0] = util::min(_min[0], (T) tri.point[2][0]);
695 _min[1] = util::min(_min[1], (T) tri.point[2][1]);
696 _min[2] = util::min(_min[2], (T) tri.point[2][2]);
697
698 _max[0] = util::max(_max[0], (T) tri.point[2][0]);
699 _max[1] = util::max(_max[1], (T) tri.point[2][1]);
700 _max[2] = util::max(_max[2], (T) tri.point[2][2]);
701 }
702
703 tri.init();
704 _triangles.push_back(tri);
705
706 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[1]),
707 _maxDist2);
708 _maxDist2 = util::max(distPoints(tri.point[2], tri.point[1]),
709 _maxDist2);
710 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[2]),
711 _maxDist2);
712 }
713}
714
715template<typename T>
716T STLmesh<T>::distPoints(Vertex<T,3>& p1, Vertex<T,3>& p2)
717{
718 return util::pow(T(p1[0] - p2[0]), 2)
719 + util::pow(T(p1[1] - p2[1]), 2)
720 + util::pow(T(p1[2] - p2[2]), 2);
721}
722
723template<typename T>
724void STLmesh<T>::print(bool full)
725{
726 if (full) {
727 int i = 0;
728 clout << "Triangles: " << std::endl;
729 typename std::vector<STLtriangle<T> >::iterator it = _triangles.begin();
730
731 for (; it != _triangles.end(); ++it) {
732 clout << i++ << ": " << it->point[0][0] << " " << it->point[0][1]
733 << " " << it->point[0][2] << " | " << it->point[1][0] << " "
734 << it->point[1][1] << " " << it->point[1][2] << " | "
735 << it->point[2][0] << " " << it->point[2][1] << " "
736 << it->point[2][2] << std::endl;
737 }
738 }
739 clout << "nTriangles=" << _triangles.size() << "; maxDist2=" << _maxDist2
740 << std::endl;
741 clout << "minPhysR(StlMesh)=(" << getMin()[0] << "," << getMin()[1] << ","
742 << getMin()[2] << ")";
743 clout << "; maxPhysR(StlMesh)=(" << getMax()[0] << "," << getMax()[1] << ","
744 << getMax()[2] << ")" << std::endl;
745}
746
747template<typename T>
748void STLmesh<T>::write(std::string fName)
749{
750 int rank = 0;
751#ifdef PARALLEL_MODE_MPI
752 rank = singleton::mpi().getRank();
753#endif
754 if (rank == 0) {
755 std::string fullName = singleton::directories().getVtkOutDir() + fName
756 + ".stl";
757 std::ofstream f(fullName.c_str());
758 f << "solid ascii " << fullName << "\n";
759
760 for (unsigned int i = 0; i < _triangles.size(); i++) {
761 f << "facet normal " << _triangles[i].normal[0] << " "
762 << _triangles[i].normal[1] << " " << _triangles[i].normal[2] << "\n";
763 f << " outer loop\n";
764 f << " vertex " << _triangles[i].point[0][0] << " "
765 << _triangles[i].point[0][1] << " " << _triangles[i].point[0][2]
766 << "\n";
767 f << " vertex " << _triangles[i].point[1][0] << " "
768 << _triangles[i].point[1][1] << " " << _triangles[i].point[1][2]
769 << "\n";
770 f << " vertex " << _triangles[i].point[2][0] << " "
771 << _triangles[i].point[2][1] << " " << _triangles[i].point[2][2]
772 << "\n";
773 f << " endloop\n";
774 f << "endfacet\n";
775 }
776 f << "endsolid\n";
777 f.close();
778 }
779 /*if (_verbose)*/clout << "Write ... OK" << std::endl;
780}
781
782template<typename T>
783bool STLmesh<T>::testRayIntersect(const std::set<unsigned int>& tris, const Vector<T,3>& pt,const Vector<T,3>& dir, Vector<T,3>& q, T& alpha)
784{
785 std::set<unsigned int>::iterator it = tris.begin();
786 for (; it != tris.end(); ++it) {
787 if (_triangles[*it].testRayIntersect(pt, dir, q, alpha) && alpha < 1) {
788 return true;
789 }
790 }
791 return false;
792}
793
794
795/*
796 * STLReader functions
797 */
798template<typename T>
799STLreader<T>::STLreader(const std::string fName, T voxelSize, T stlSize,
800 RayMode method, bool verbose, T overlap, T max)
801 : _voxelSize(voxelSize),
802 _stlSize(stlSize),
803 _overlap(overlap),
804 _fName(fName),
805 _mesh(fName, stlSize),
806 _verbose(verbose),
807 clout(std::cout, "STLreader")
808{
809 this->getName() = "STLreader";
810
811 if (_verbose) {
812 clout << "Voxelizing ..." << std::endl;
813 }
814
815 Vector<T,3> extension = _mesh.getMax() - _mesh.getMin();
816 if ( util::nearZero(max) ) {
817 max = util::max(extension[0], util::max(extension[1], extension[2])) + _voxelSize;
818 }
819 int j = 0;
820 for (; _voxelSize * util::pow(2, j) < max; j++)
821 ;
822 Vector<T,3> center;
823 T radius = _voxelSize * util::pow(2, j - 1);
824
826 for (unsigned i = 0; i < 3; i++) {
827 center[i] = (_mesh.getMin()[i] + _mesh.getMax()[i]) / 2. - _voxelSize / 4.;
828 }
829
831 _tree = new Octree<T>(center, radius, &_mesh, j, _overlap);
832
834 for (int i = 0; i < 3; i++) {
835 this->_myMin[i] = center[i] + _voxelSize / 2.;
836 this->_myMax[i] = center[i] - _voxelSize / 2.;
837 }
838 for (int i = 0; i < 3; i++) {
839 while (this->_myMin[i] > _mesh.getMin()[i]) {
840 this->_myMin[i] -= _voxelSize;
841 }
842 while (this->_myMax[i] < _mesh.getMax()[i]) {
843 this->_myMax[i] += _voxelSize;
844 }
845 this->_myMax[i] -= _voxelSize;
846 this->_myMin[i] += _voxelSize;
847 }
848
850 switch (method) {
851 case RayMode::Robust:
852 indicate1();
853 break;
854 case RayMode::DoubleRay:
855 indicate3();
856 break;
857 case RayMode::FastRayX:
858 indicate2_Xray();
859 break;
860 case RayMode::FastRayY:
861 indicate2_Yray();
862 break;
863 default:
864 indicate2();
865 break;
866 }
867
868 if (_verbose) {
869 print();
870 }
871 if (_verbose) {
872 clout << "Voxelizing ... OK" << std::endl;
873 }
874}
875
876/*
877 * STLReader functions
878 */
879template<typename T>
880STLreader<T>::STLreader(const std::vector<std::vector<T>> meshPoints, T voxelSize, T stlSize,
881 RayMode method, bool verbose, T overlap, T max)
882 : _voxelSize(voxelSize),
883 _stlSize(stlSize),
884 _overlap(overlap),
885 _fName("meshPoints.stl"),
886 _mesh(meshPoints, stlSize),
887 _verbose(verbose),
888 clout(std::cout, "STLreader")
889{
890 this->getName() = "STLreader";
891
892 if (_verbose) {
893 clout << "Voxelizing ..." << std::endl;
894 }
895
896 Vector<T,3> extension = _mesh.getMax() - _mesh.getMin();
897 if ( util::nearZero(max) ) {
898 max = util::max(extension[0], util::max(extension[1], extension[2])) + _voxelSize;
899 }
900 int j = 0;
901 for (; _voxelSize * util::pow(2, j) < max; j++)
902 ;
903 Vector<T,3> center;
904 T radius = _voxelSize * util::pow(2, j - 1);
905
907 for (unsigned i = 0; i < 3; i++) {
908 center[i] = (_mesh.getMin()[i] + _mesh.getMax()[i]) / 2. - _voxelSize / 4.;
909 }
910
912
913 _tree = new Octree<T>(center, radius, &_mesh, j, _overlap);
914
916 for (int i = 0; i < 3; i++) {
917 this->_myMin[i] = center[i] + _voxelSize / 2.;
918 this->_myMax[i] = center[i] - _voxelSize / 2.;
919 }
920 for (int i = 0; i < 3; i++) {
921 while (this->_myMin[i] > _mesh.getMin()[i]) {
922 this->_myMin[i] -= _voxelSize;
923 }
924 while (this->_myMax[i] < _mesh.getMax()[i]) {
925 this->_myMax[i] += _voxelSize;
926 }
927 this->_myMax[i] -= _voxelSize;
928 this->_myMin[i] += _voxelSize;
929 }
930 //automaticly choose the method with minimum extension in its direction
931
932 /*if(extension[0] == std::min_element(extension.begin(), extension.end())){
933 method = 4;
934 }
935 else if(extension[1] == std::min_element(extension.begin(), extension.end())){
936 method = 5;
937 }
938 else if(extension[2] == std::min_element(extension.begin(), extension.end())){
939 method = 0;
940 }
941 */
942
943
944 // Indicate nodes of the tree. (Inside/Outside)
945 switch (method) {
946 case RayMode::Robust:
947 indicate1();
948 break;
949 case RayMode::DoubleRay:
950 indicate3();
951 break;
952 case RayMode::FastRayX:
953 indicate2_Xray();
954 break;
955 case RayMode::FastRayY:
956 indicate2_Yray();
957 break;
958 default:
959 indicate2();
960 break;
961 }
962
963 setNormalsOutside();
964
965 if (_verbose) {
966 print();
967 }
968 if (_verbose) {
969 clout << "Voxelizing ... OK" << std::endl;
970 }
971}
972
973template<typename T>
974STLreader<T>::~STLreader()
975{
976 delete _tree;
977}
978
979/*
980 * Old indicate function (slower, more stable)
981 * Define three rays (X-, Y-, Z-direction) for each leaf and count intersections
982 * with STL for each ray. Odd number of intersection means inside (Majority vote).
983 */
984
985template<typename T>
986void STLreader<T>::indicate1()
987{
988 std::vector<Octree<T>*> leafs;
989 _tree->getLeafs(leafs);
990 typename std::vector<Octree<T>*>::iterator it = leafs.begin();
991 Vector<T,3> dir, pt, s;
992
993 int intersections = 0;
994 int inside = 0;
995 Octree<T>* node = nullptr;
996 T step = 1. / 1000. * _voxelSize;
997 for (; it != leafs.end(); ++it) {
998 inside = 0;
999
1000 pt = (*it)->getCenter();
1001 intersections = 0;
1002 s = pt; // + step;
1003
1005 dir[0] = 1;
1006 dir[1] = 0;
1007 dir[2] = 0;
1008 while (s[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1009 node = _tree->find(s, (*it)->getMaxdepth());
1010 intersections += node->testIntersection(pt, dir);
1011 node->intersectRayNode(pt, dir, s);
1012 s = s + step * dir;
1013 }
1014 inside += (intersections % 2);
1015
1017 intersections = 0;
1018 s = pt; // + step;
1019 dir[0] = 0;
1020 dir[1] = 1;
1021 dir[2] = 0;
1022 while (s[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1023 node = _tree->find(s, (*it)->getMaxdepth());
1024 intersections += node->testIntersection(pt, dir);
1025 node->intersectRayNode(pt, dir, s);
1026 s = s + step * dir;
1027 }
1028 inside += (intersections % 2);
1029
1031 intersections = 0;
1032 s = pt; // + step;
1033 dir[0] = 0;
1034 dir[1] = 0;
1035 dir[2] = 1;
1036 while (s[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1037 node = _tree->find(s, (*it)->getMaxdepth());
1038 intersections += node->testIntersection(pt, dir);
1039 node->intersectRayNode(pt, dir, s);
1040 s = s + step * dir;
1041 }
1042 inside += (intersections % 2);
1043 (*it)->setInside(inside > 1);
1044 }
1045}
1046
1047/*
1048 * New indicate function (faster, less stable)
1049 * Define ray in Z-direction for each Voxel in XY-layer. Indicate all nodes on the fly.
1050 */
1051template<typename T>
1052void STLreader<T>::indicate2()
1053{
1054 T rad = _tree->getRadius();
1055 Vector<T,3> rayPt = _tree->getCenter() - rad + .5 * _voxelSize;
1056 Vector<T,3> pt = rayPt;
1057 Vector<T,3> rayDir;
1058 rayDir[0] = 0.;
1059 rayDir[1] = 0.;
1060 rayDir[2] = 1.;
1061 //Vector<T,3> maxEdge = _tree->getCenter() + rad;
1062
1063 T step = 1. / 1000. * _voxelSize;
1064
1065 Octree<T>* node = nullptr;
1066 unsigned short rayInside = 0;
1067 Vector<T,3> nodeInters;
1068 while (pt[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1069 node = _tree->find(pt);
1070 nodeInters = pt;
1071 nodeInters[2] = node->getCenter()[2] - node->getRadius();
1072 rayInside = 0;
1073 while (pt[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1074 node = _tree->find(pt);
1075 nodeInters = pt;
1076 nodeInters[2] = node->getCenter()[2] - node->getRadius();
1077 rayInside = 0;
1078 while (pt[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1079 node = _tree->find(pt);
1080 node->checkRay(nodeInters, rayDir, rayInside);
1081 node->intersectRayNode(pt, rayDir, nodeInters);
1082 pt = nodeInters + step * rayDir;
1083 }
1084 pt[2] = rayPt[2];
1085 pt[1] += _voxelSize;
1086 }
1087 pt[1] = rayPt[1];
1088 pt[0] += _voxelSize;
1089 }
1090}
1091
1092
1093/*
1094 * New indicate function (faster, less stable)
1095 * Define ray in X-direction for each Voxel in YZ-layer. Indicate all nodes on the fly.
1096 */
1097
1098template<typename T>
1099void STLreader<T>::indicate2_Xray()
1100{
1101 T rad = _tree->getRadius();
1102 Vector<T,3> rayPt = _tree->getCenter() - rad + .5 * _voxelSize;
1103 Vector<T,3> pt = rayPt;
1104 Vector<T,3> rayDir;
1105 rayDir[0] = 1.;
1106 rayDir[1] = 0.;
1107 rayDir[2] = 0.;
1108 //Vector<T,3> maxEdge = _tree->getCenter() + rad;
1109
1110 T step = 1. / 1000. * _voxelSize;
1111
1112 Octree<T>* node = nullptr;
1113 unsigned short rayInside = 0;
1114 Vector<T,3> nodeInters;
1115 while (pt[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1116 node = _tree->find(pt);
1117 nodeInters = pt;
1118 nodeInters[0] = node->getCenter()[0] - node->getRadius();
1119 rayInside = 0;
1120 while (pt[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1121 node = _tree->find(pt);
1122 nodeInters = pt;
1123 nodeInters[0] = node->getCenter()[0] - node->getRadius();
1124 rayInside = 0;
1125 while (pt[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1126 node = _tree->find(pt);
1127 node->checkRay(nodeInters, rayDir, rayInside);
1128 node->intersectRayNode(pt, rayDir, nodeInters);
1129 pt = nodeInters + step * rayDir;
1130 }
1131 pt[0] = rayPt[0];
1132 pt[1] += _voxelSize;
1133 }
1134 pt[1] = rayPt[1];
1135 pt[2] += _voxelSize;
1136 }
1137}
1138
1139/*
1140 * New indicate function (faster, less stable)
1141 * Define ray in Y-direction for each Voxel in XZ-layer. Indicate all nodes on the fly.
1142 */
1143
1144template<typename T>
1145void STLreader<T>::indicate2_Yray()
1146{
1147 T rad = _tree->getRadius();
1148 Vector<T,3> rayPt = _tree->getCenter() - rad + .5 * _voxelSize;
1149 Vector<T,3> pt = rayPt;
1150 Vector<T,3> rayDir;
1151 rayDir[0] = 0.;
1152 rayDir[1] = 1.;
1153 rayDir[2] = 0.;
1154 //Vector<T,3> maxEdge = _tree->getCenter() + rad;
1155
1156 T step = 1. / 1000. * _voxelSize;
1157
1158 Octree<T>* node = nullptr;
1159 unsigned short rayInside = 0;
1160 Vector<T,3> nodeInters;
1161 while (pt[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1162 node = _tree->find(pt);
1163 nodeInters = pt;
1164 nodeInters[1] = node->getCenter()[1] - node->getRadius();
1165 rayInside = 0;
1166 while (pt[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1167 node = _tree->find(pt);
1168 nodeInters = pt;
1169 nodeInters[1] = node->getCenter()[1] - node->getRadius();
1170 rayInside = 0;
1171 while (pt[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1172 node = _tree->find(pt);
1173 node->checkRay(nodeInters, rayDir, rayInside);
1174 node->intersectRayNode(pt, rayDir, nodeInters);
1175 pt = nodeInters + step * rayDir;
1176 }
1177 pt[1] = rayPt[1];
1178 pt[0] += _voxelSize;
1179 }
1180 pt[0] = rayPt[0];
1181 pt[2] += _voxelSize;
1182 }
1183}
1184
1185/*
1186 * Double ray approach: two times (X-, Y-, Z-direction) for each leaf.
1187 * Could be use to deal with double layer triangles and face intersections.
1188 */
1189template<typename T>
1190void STLreader<T>::indicate3()
1191{
1192 std::vector<Octree<T>*> leafs;
1193 _tree->getLeafs(leafs);
1194 typename std::vector<Octree<T>*>::iterator it = leafs.begin();
1195
1196 Vector<T,3> dir, pt, s;
1197 Octree<T>* node = nullptr;
1198 T step = 1. / 1000. * _voxelSize;
1199 int intersections;
1200 int sum_intersections;
1201
1202 for (; it != leafs.end(); ++it) {
1203 pt = (*it)->getCenter();
1204 intersections = 0;
1205 sum_intersections = 0;
1206 s = pt; // + step;
1207
1209 dir[0] = 1;
1210 dir[1] = 0;
1211 dir[2] = 0;
1212 while (s[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1213 node = _tree->find(s, (*it)->getMaxdepth());
1214 intersections = node->testIntersection(pt, dir);
1215 node->intersectRayNode(pt, dir, s);
1216 s = s + step * dir;
1217 if (intersections > 0) {
1218 sum_intersections++;
1219 break;
1220 }
1221 }
1222
1224 intersections = 0;
1225 s = pt; // + step;
1226 dir[0] = 0;
1227 dir[1] = 1;
1228 dir[2] = 0;
1229 while (s[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1230 node = _tree->find(s, (*it)->getMaxdepth());
1231 intersections = node->testIntersection(pt, dir);
1232 node->intersectRayNode(pt, dir, s);
1233 s = s + step * dir;
1234 if (intersections > 0) {
1235 sum_intersections++;
1236 break;
1237 }
1238 }
1239
1241 intersections = 0;
1242 s = pt; // + step;
1243 dir[0] = 0;
1244 dir[1] = 0;
1245 dir[2] = 1;
1246 while (s[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1247 node = _tree->find(s, (*it)->getMaxdepth());
1248 intersections = node->testIntersection(pt, dir);
1249 node->intersectRayNode(pt, dir, s);
1250 s = s + step * dir;
1251 if (intersections > 0) {
1252 sum_intersections++;
1253 break;
1254 }
1255 }
1256
1258 intersections = 0;
1259 s = pt; // + step;
1260 dir[0] = -1;
1261 dir[1] = 0;
1262 dir[2] = 0;
1263 while (s[0] > _mesh.getMin()[0] - std::numeric_limits<T>::epsilon()) {
1264 node = _tree->find(s, (*it)->getMaxdepth());
1265 intersections = node->testIntersection(pt, dir);
1266 node->intersectRayNode(pt, dir, s);
1267 s = s + step * dir;
1268 if (intersections > 0) {
1269 sum_intersections++;
1270 break;
1271 }
1272 }
1273
1275 intersections = 0;
1276 s = pt; // + step;
1277 dir[0] = 0;
1278 dir[1] = -1;
1279 dir[2] = 0;
1280 while (s[1] > _mesh.getMin()[1] - std::numeric_limits<T>::epsilon()) {
1281 node = _tree->find(s, (*it)->getMaxdepth());
1282 intersections = node->testIntersection(pt, dir);
1283 node->intersectRayNode(pt, dir, s);
1284 s = s + step * dir;
1285 if (intersections > 0) {
1286 sum_intersections++;
1287 break;
1288 }
1289 }
1290
1292 intersections = 0;
1293 s = pt; // + step;
1294 dir[0] = 0;
1295 dir[1] = 0;
1296 dir[2] = -1;
1297 while (s[2] > _mesh.getMin()[2] - std::numeric_limits<T>::epsilon()) {
1298 node = _tree->find(s, (*it)->getMaxdepth());
1299 intersections = node->testIntersection(pt, dir);
1300 node->intersectRayNode(pt, dir, s);
1301 s = s + step * dir;
1302 if (intersections > 0) {
1303 sum_intersections++;
1304 break;
1305 }
1306 }
1307 (*it)->setInside(sum_intersections > 5);
1308 }
1309}
1310
1311template<typename T>
1312bool STLreader<T>::operator() (bool output[], const T input[])
1313{
1314 output[0] = false;
1315 if (isInsideRootTree(input)) {
1316 std::vector<T> tmp(input, input + 3);
1317 output[0] = _tree->find(tmp)->getInside();
1318 }
1319 return true;
1320}
1321
1322
1323template<typename T>
1324bool STLreader<T>::isInsideRootTree(const T input[])
1325{
1326 T coords = _tree->getRadius();
1327 Vector<T,3> c(_tree->getCenter());
1328 return c[0] - coords < input[0] && input[0] < c[0] + coords && c[1] - coords < input[1]
1329 && input[1] < c[1] + coords && c[2] - coords < input[2] && input[2] < c[2] + coords;
1330}
1331
1332
1333template<typename T>
1334Vector<T,3> STLreader<T>::closestPointInBoundingBox(const Vector<T,3>& input)
1335{
1336 T coords = _tree->getRadius();
1337 Vector<T,3> c(_tree->getCenter());
1338 Vector<T,3> closestPoint;
1339 for(int i = 0; i < 3; ++i) {
1340 closestPoint[i] = util::max(c[i] - coords, util::min(c[i] + coords, input[i]));
1341 }
1342 return closestPoint;
1343}
1344
1345
1346template<typename T>
1347bool STLreader<T>::distance(T& distance, const Vector<T,3>& origin,
1348 const Vector<T,3>& direction, int iC)
1349{
1350 Octree<T>* node = nullptr;
1351 Vector<T,3> dir(direction);
1352 dir = normalize(dir);
1353 Vector<T,3> extends = _mesh.getMax() - _mesh.getMin();
1354 Vector<T,3> pt(origin);
1355 Vector<T,3> q;
1356 Vector<T,3> s;
1357 Vector<T,3> center = _mesh.getMin() + 1 / 2. * extends;
1358 T step = _voxelSize / 1000., a = 0;
1359
1360 for (int i = 0; i < 3; i++) {
1361 extends[i] /= 2.;
1362 }
1363
1364 if (!(_mesh.getMin()[0] < origin[0] && origin[0] < _mesh.getMax()[0]
1365 && _mesh.getMin()[1] < origin[1] && origin[1] < _mesh.getMax()[1]
1366 && _mesh.getMin()[2] < origin[2] && origin[2] < _mesh.getMax()[2])) {
1367 T t = T(), d = T();
1368 bool foundQ = false;
1369
1370 if (dir[0] > 0) {
1371 d = _mesh.getMin()[0];
1372 t = (d - origin[0]) / dir[0];
1373 pt[0] = origin[0] + (t + step) * dir[0];
1374 pt[1] = origin[1] + (t + step) * dir[1];
1375 pt[2] = origin[2] + (t + step) * dir[2];
1376
1377 if (_mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]
1378 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1379 foundQ = true;
1380 }
1381 }
1382 else if (dir[0] < 0) {
1383 d = _mesh.getMax()[0];
1384 t = (d - origin[0]) / dir[0];
1385 pt[0] = origin[0] + (t + step) * dir[0];
1386 pt[1] = origin[1] + (t + step) * dir[1];
1387 pt[2] = origin[2] + (t + step) * dir[2];
1388 if (_mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]
1389 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1390 foundQ = true;
1391 }
1392 }
1393
1394 if (dir[1] > 0 && !foundQ) {
1395 d = _mesh.getMin()[1];
1396 t = (d - origin[1]) / dir[1];
1397 pt[0] = origin[0] + (t + step) * dir[0];
1398 pt[1] = origin[1] + (t + step) * dir[1];
1399 pt[2] = origin[2] + (t + step) * dir[2];
1400 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1401 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1402 foundQ = true;
1403 }
1404 }
1405 else if (dir[1] < 0 && !foundQ) {
1406 d = _mesh.getMax()[1];
1407 t = (d - origin[1]) / dir[1];
1408 pt[0] = origin[0] + (t + step) * dir[0];
1409 pt[1] = origin[1] + (t + step) * dir[1];
1410 pt[2] = origin[2] + (t + step) * dir[2];
1411 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1412 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1413 foundQ = true;
1414 }
1415 }
1416
1417 if (dir[2] > 0 && !foundQ) {
1418 d = _mesh.getMin()[2];
1419 t = (d - origin[2]) / dir[2];
1420 pt[0] = origin[0] + (t + step) * dir[0];
1421 pt[1] = origin[1] + (t + step) * dir[1];
1422 pt[2] = origin[2] + (t + step) * dir[2];
1423 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1424 && _mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]) {
1425 foundQ = true;
1426 }
1427 }
1428 else if (dir[2] < 0 && !foundQ) {
1429 d = _mesh.getMax()[2];
1430 t = (d - origin[2]) / dir[2];
1431 pt[0] = origin[0] + (t + step) * dir[0];
1432 pt[1] = origin[1] + (t + step) * dir[1];
1433 pt[2] = origin[2] + (t + step) * dir[2];
1434 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1435 && _mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]) {
1436 foundQ = true;
1437 }
1438 }
1439
1440 if (!foundQ) {
1441 return false;
1442 }
1443 }
1444
1445 while ((util::fabs(pt[0] - center[0]) < extends[0])
1446 && (util::fabs(pt[1] - center[1]) < extends[1])
1447 && (util::fabs(pt[2] - center[2]) < extends[2])) {
1448 node = _tree->find(pt);
1449 if (node->closestIntersection(Vector<T,3>(origin), dir, q, a)) {
1450 Vector<T,3> vek(q - Vector<T,3>(origin));
1451 distance = norm(vek);
1452 return true;
1453 }
1454 else {
1455 Octree<T>* tmpNode = _tree->find(pt);
1456 tmpNode->intersectRayNode(pt, dir, s);
1457 for (int i = 0; i < 3; i++) {
1458 pt[i] = s[i] + step * dir[i];
1459 }
1460 }
1461 }
1462
1463 return false;
1464}
1465
1466template<typename T>
1467template<typename F>
1468void STLreader<T>::iterateOverCloseTriangles(const PhysR<T,3>& pt, F func, Octree<T>* leafNode) {
1469 // Find the leaf node in the Octree that contains the input point
1470 leafNode = _tree->find(pt);
1471
1472 if (leafNode && !leafNode->getTriangles().empty()) {
1473 const std::vector<unsigned int>& triangleIndices = leafNode->getTriangles();
1474 for (unsigned int idx : triangleIndices) {
1475 const STLtriangle<T>& triangle = _mesh.getTri(idx);
1476 func(triangle);
1477 }
1478 }
1479 else {
1480 for (const STLtriangle<T>& triangle : _mesh.getTriangles()) {
1481 func(triangle);
1482 }
1483 }
1484};
1485
1486template<typename T>
1487Vector<T,3> STLreader<T>::evalNormalOnSurface(const PhysR<T,3>& pt, const Vector<T,3>& fallbackNormal)
1488{
1489 // Check if the position is on the corner of a triangle
1490 unsigned countTriangles = 0;
1491 Vector<T,3> normal(T(0));
1492
1493 iterateOverCloseTriangles(pt, [&](const STLtriangle<T>& triangle){
1494 // TODO: Calculate angle-weighted psuedonormal (see 10.1109/TVCG.2005.49) in case the point lies on corners
1495 // Edges correspond to an unweighted average (as calculated below) anyway
1496 if (triangle.isPointInside(pt)) {
1497 ++countTriangles;
1498 normal+=triangle.getNormal();
1499 }
1500 });
1501
1502 if (countTriangles > 0) {
1503 return normal / countTriangles;
1504 }
1505
1506 // If the provided point isn't located on the surface, return the predefined fallback
1507 return fallbackNormal;
1508}
1509
1510template<typename T>
1511template<SignMode SIGNMODE>
1512Vector<T,3> STLreader<T>::evalSurfaceNormal(const Vector<T,3>& origin)
1513{
1514 Vector<T,3> normal(0.);
1515 Vector<T,3> closestPointOnSurface(0.);
1516 const STLtriangle<T>* closestTriangle = nullptr;
1517 T distance = std::numeric_limits<T>::max();
1518 const PhysR<T,3> pt(closestPointInBoundingBox(origin));
1519
1520 iterateOverCloseTriangles(pt, [&](const STLtriangle<T>& triangle) {
1521 PhysR<T,3> const pointOnTriangle = triangle.closestPtPointTriangle(pt);
1522 PhysR<T,3> const currDistance = pt - pointOnTriangle;
1523 T currDistanceNorm = norm_squared(currDistance);
1524 if (distance > currDistanceNorm) {
1525 normal = currDistance;
1526 distance = currDistanceNorm;
1527 closestPointOnSurface = pointOnTriangle;
1528 closestTriangle = &triangle;
1529 }
1530 });
1531
1532 distance = util::sqrt(distance);
1533 if (!util::nearZero(distance)) {
1535 normal *= evalSignForSignedDistance<SIGNMODE>(origin, distance);
1536 }
1537 else {
1538 normal = evalNormalOnSurface(closestPointOnSurface, closestTriangle->getNormal());
1539 // The below isn't necessary because all normals are set to point outside by default (see constructor)
1540
1541 /* const T tmpStepSize = _voxelSize * T{0.01}; */
1542 /* const short signBefore = evalSignForSignedDistance(origin - tmpStepSize * normal, *closestTriangle); */
1543 /* const short signAfter = evalSignForSignedDistance(origin + tmpStepSize * normal, *closestTriangle); */
1544
1545 /* if(signBefore > 0 || signAfter < 0) { */
1546 /* normal *= T{-1}; */
1547 /* } */
1548 }
1549
1550 return normal;
1551}
1552
1553
1554template<typename T>
1555template<SignMode SIGNMODE>
1556Vector<T,3> STLreader<T>::evalSurfaceNormalForPseudoNormal(const Vector<T,3>& origin, Vector<T,3> & outputPointOnSurface)
1557{
1558 Vector<T,3> normal(0.);
1559 Vector<T,3> closestPointOnSurface(0.);
1560 const STLtriangle<T>* closestTriangle = nullptr;
1561 T distance = std::numeric_limits<T>::max();
1562 const PhysR<T,3> pt(closestPointInBoundingBox(origin));
1563
1564 iterateOverCloseTriangles(pt, [&](const STLtriangle<T>& triangle) {
1565 PhysR<T,3> const pointOnTriangle = triangle.closestPtPointTriangle(pt);
1566 PhysR<T,3> const currDistance = pt - pointOnTriangle;
1567 T currDistanceNorm = norm_squared(currDistance);
1568 if (distance > currDistanceNorm) {
1569 normal = currDistance;
1570 distance = currDistanceNorm;
1571 closestPointOnSurface = pointOnTriangle;
1572 closestTriangle = &triangle;
1573 }
1574 });
1575
1576/* distance = util::sqrt(distance);
1577 if (!util::nearZero(distance)) {
1578 normal = normal/distance;
1579 normal *= evalSignForSignedDistance<SIGNMODE>(origin, distance);
1580 }
1581 else { */
1582 normal = evalNormalOnSurface(closestPointOnSurface, closestTriangle->getNormal());
1583 outputPointOnSurface = closestPointOnSurface;
1584
1585 // The below isn't necessary because all normals are set to point outside by default (see constructor)
1586
1587 /* const T tmpStepSize = _voxelSize * T{0.01}; */
1588 /* const short signBefore = evalSignForSignedDistance(origin - tmpStepSize * normal, *closestTriangle); */
1589 /* const short signAfter = evalSignForSignedDistance(origin + tmpStepSize * normal, *closestTriangle); */
1590
1591 /* if(signBefore > 0 || signAfter < 0) { */
1592 /* normal *= T{-1}; */
1593 /* } */
1594 //}
1595
1596 return normal;
1597}
1598
1599
1600template <typename T>
1601template <SignMode SIGNMODE>
1602short STLreader<T>::evalSignForSignedDistance(
1603 const Vector<T,3>& pt, [[maybe_unused]] const T distance, Vector<T,3> vecdist, STLtriangle<T> stlT)
1604{
1605 short sign;
1606
1607 if constexpr (SIGNMODE == SignMode::EXACT)
1608 {
1609
1610 if(distance < _voxelSize) {
1611
1612 sign = evalSignForSignedDistanceFromWindingNumber(pt);
1613 //sign = evalSignForSignedDistanceFromNormal(stlT.getNormal(), vecdist);
1614 /* Vector<T,3> closestPointOnSurface (0.);
1615 Vector<T,3> normal = evalSurfaceNormalForPseudoNormal<SignMode::EXACT>(pt,closestPointOnSurface);
1616 sign = evalSignForSignedDistanceFromNormal(normal, pt - closestPointOnSurface);*/
1617 }
1618 else
1619 {
1620 sign = evalSignForSignedDistanceFromCache(pt);
1621 }
1622 }
1623 else {
1624 sign = evalSignForSignedDistanceFromCache(pt);
1625 }
1626 return sign;
1627}
1628
1629template<typename T>
1630short STLreader<T>::evalSignForSignedDistanceFromPseudonormal(const Vector<T,3>& pseudonormal, const Vector<T,3>& distance)
1631{
1632 const T projection = pseudonormal * distance;
1633
1634 if(projection > 0) {
1635 return 1;
1636 }
1637 else if (projection < 0 ) {
1638 return -1;
1639 }
1640 return 0;
1641}
1642
1643
1644//aproximation not working for concave surfaces
1645template<typename T>
1646short STLreader<T>::evalSignForSignedDistanceFromNormal(const Vector<T,3>& normal, const Vector<T,3>& distance)
1647{
1648 const T projection = normal * distance;
1649
1650 if(projection > 0) {
1651 return 1;
1652 }
1653 else if (projection < 0 ) {
1654 return -1;
1655 }
1656 return 0;
1657}
1658
1659
1660
1661
1662template<typename T>
1663short STLreader<T>::evalSignForSignedDistanceFromWindingNumber(const Vector<T,3>& pt)
1664{
1665 T windingNumber{0};
1666
1667 for (const STLtriangle<T>& triangle : _mesh.getTriangles()) {
1668 const PhysR<T,3> a = triangle.point[0] - pt;
1669 const PhysR<T,3> b = triangle.point[1] - pt;
1670 const PhysR<T,3> c = triangle.point[2] - pt;
1671
1672 const T aNorm = norm(a);
1673 const T bNorm = norm(b);
1674 const T cNorm = norm(c);
1675
1676 const T numerator = a * crossProduct3D(b, c);
1677 const T denominator = aNorm * bNorm * cNorm + (a*b) * cNorm + (b*c) * aNorm + (c*a) * bNorm;
1678
1679 windingNumber += util::atan2(numerator,denominator);
1680 }
1681 windingNumber *= 2;
1682
1683 if(int(util::round(windingNumber))%2 == 0) {
1684 return 1;
1685 }
1686 return -1;
1687}
1688
1689
1690template <typename T>
1691T STLreader<T>::signedDistance(const Vector<T, 3>& input)
1692{
1693 return signedDistance<SignMode::CACHED>(input);
1694}
1695
1696template <typename T>
1697T STLreader<T>::signedDistanceExact(const Vector<T, 3>& input)
1698{
1699 return signedDistance<SignMode::EXACT>(input);
1700}
1701
1702
1703
1704
1705template <typename T>
1706template <SignMode SIGNMODE>
1707T STLreader<T>::signedDistance(const Vector<T, 3>& input)
1708{
1709 T distanceNorm = std::numeric_limits<T>::max();
1710 //Vector<T,3> distVec (std::numeric_limits<T>::max(),std::numeric_limits<T>::max(),std::numeric_limits<T>::max());
1711 T extraDistance{};
1712 //STLtriangle<T> stlT;
1713
1714 const auto evalDistance = [&](const Vector<T,3>& pt) {
1715 iterateOverCloseTriangles(pt, [&](const STLtriangle<T>& triangle) {
1716 const PhysR<T, 3> currPointOnTriangle =
1717 triangle.closestPtPointTriangle(pt);
1718 const Vector<T, 3> currDistance = pt - currPointOnTriangle;
1719 T currDistanceNorm = norm_squared(currDistance);
1720 distanceNorm = util::min(distanceNorm, currDistanceNorm);
1721 /* if(distanceNorm -currDistanceNorm < std::numeric_limits<T>::min())
1722 {distVec = currDistanceNorm;
1723 stlT = triangle;}*/
1724 });
1725 return util::sqrt(distanceNorm);
1726 };
1727
1728 if(!isInsideRootTree(input.data())) {
1729 const PhysR<T,3> ptInsideTree = closestPointInBoundingBox(input);
1730 const T distance = evalDistance(ptInsideTree);
1731 extraDistance = norm(ptInsideTree - input);
1732 return distance + extraDistance;
1733 }
1734 const T distance = evalDistance(input);
1735
1736 return evalSignForSignedDistance<SIGNMODE>(input, distance) * distance;
1737}
1738
1739
1740template<typename T>
1741short STLreader<T>::evalSignForSignedDistanceFromCache(const Vector<T,3>& pt)
1742{
1743 bool isInside;
1744 this->operator()(&isInside, pt.data());
1745 return (isInside ? -1 : 1);
1746}
1747
1748
1749
1750/* template <typename T>
1751template <SignMode SIGNMODE>
1752T STLreader<T>::signedDistance(const Vector<T, 3>& input)
1753{
1754 T distanceNorm = std::numeric_limits<T>::max();
1755 T extraDistance{};
1756
1757 const auto evalDistance = [&](const Vector<T,3>& pt) {
1758 iterateOverCloseTriangles(pt, [&](const STLtriangle<T>& triangle) {
1759 const PhysR<T, 3> currPointOnTriangle =
1760 triangle.closestPtPointTriangle(pt);
1761 const Vector<T, 3> currDistance = pt - currPointOnTriangle;
1762 T currDistanceNorm = norm_squared(currDistance);
1763 distanceNorm = util::min(distanceNorm, currDistanceNorm);
1764 });
1765 return util::sqrt(distanceNorm);
1766 };
1767
1768 if(!isInsideRootTree(input.data())) {
1769 const PhysR<T,3> ptInsideTree = closestPointInBoundingBox(input);
1770 const T distance = evalDistance(ptInsideTree);
1771 extraDistance = norm(ptInsideTree - input);
1772 return distance + extraDistance;
1773 }
1774 const T distance = evalDistance(input);
1775
1776 return evalSignForSignedDistance<SIGNMODE>(input, distance) * distance;
1777}*/
1778
1779template <typename T>
1780Vector<T,3> STLreader<T>::surfaceNormal(const Vector<T,3>& pos, const T meshSize)
1781{
1782 return surfaceNormal<SignMode::CACHED>(pos, meshSize);
1783}
1784
1785template <typename T>
1786Vector<T,3> STLreader<T>::surfaceNormalExact(const Vector<T,3>& pos, const T meshSize)
1787{
1788 return surfaceNormal<SignMode::EXACT>(pos, meshSize);
1789}
1790
1791template <typename T>
1792template <SignMode SIGNMODE>
1793Vector<T,3> STLreader<T>::surfaceNormal(const Vector<T,3>& pos, const T meshSize)
1794{
1795 return evalSurfaceNormal<SIGNMODE>(pos);
1796}
1797
1798template<typename T>
1799void STLreader<T>::print()
1800{
1801 _mesh.print();
1802 _tree->print();
1803 clout << "voxelSize=" << _voxelSize << "; stlSize=" << _stlSize << std::endl;
1804 clout << "minPhysR(VoxelMesh)=(" << this->_myMin[0] << "," << this->_myMin[1]
1805 << "," << this->_myMin[2] << ")";
1806 clout << "; maxPhysR(VoxelMesh)=(" << this->_myMax[0] << ","
1807 << this->_myMax[1] << "," << this->_myMax[2] << ")" << std::endl;
1808}
1809
1810template<typename T>
1811void STLreader<T>::writeOctree()
1812{
1813 _tree->write(_fName);
1814}
1815
1816template<typename T>
1817void STLreader<T>::writeSTL(std::string stlName)
1818{
1819 if (stlName == "") {
1820 _mesh.write(_fName);
1821 }
1822 else {
1823 _mesh.write(stlName);
1824 }
1825}
1826
1827template<typename T>
1828void STLreader<T>::setNormalsOutside()
1829{
1830 unsigned int noTris = _mesh.triangleSize();
1831 Vector<T,3> center;
1832 //Octree<T>* node = nullptr;
1833 for (unsigned int i = 0; i < noTris; i++) {
1834 center = _mesh.getTri(i).getCenter();
1835 if (_tree->find(
1836 center + _mesh.getTri(i).normal * util::sqrt(3.) * _voxelSize)->getInside()) {
1837 // cout << "Wrong direction" << std::endl;
1838 Vector<T,3> pt(_mesh.getTri(i).point[0]);
1839 _mesh.getTri(i).point[0] = _mesh.getTri(i).point[2];
1840 _mesh.getTri(i).point[2] = pt;
1841 _mesh.getTri(i).init();
1842 // _mesh.getTri(i).getNormal()[0] *= -1.;
1843 // _mesh.getTri(i).getNormal()[1] *= -1.;
1844 // _mesh.getTri(i).getNormal()[2] *= -1.;
1845 }
1846 }
1847}
1848
1849template<typename T>
1850void STLreader<T>::setBoundaryInsideNodes()
1851{
1852 std::vector<Octree<T>*> leafs;
1853 _tree->getLeafs(leafs);
1854 for (auto it = leafs.begin(); it != leafs.end(); ++it) {
1855 if ((*it)->getBoundaryNode()) {
1856 (*it)->setInside(true);
1857 }
1858 }
1859}
1860
1861} // namespace olb
1862
1863#endif
class for marking output with some text
STLmesh(std::string, T stlSize=1.)
Constructs a new STLmesh from a file.
Definition stlReader.hh:452
Wrapper functions that simplify the use of MPI.
platform_constant int c[Q][D]
Definition functions.h:57
platform_constant Fraction s[Q]
Definition mrt.h:65
platform_constant Fraction t[Q]
Definition rtlbm.h:45
constexpr int q() any_platform
Definition functions.h:140
constexpr int d() any_platform
Definition functions.h:127
T triangle(Vector< T, 2 > p, Vector< T, 2 > a, Vector< T, 2 > b, Vector< T, 2 > c) any_platform
Exact signed distance to the surface of two-dimensional triangle.
Definition sdf.h:152
Distribution< T > normal(T mean, T stddev)
bool distance(S &distance, const Vector< S, D > &origin, const Vector< S, D > &direction, S precision, S pitch, F1 isInside, F2 isInsideBoundingBox)
Expr sqrt(Expr x)
Definition expr.cpp:225
int sign(T val) any_platform
Definition util.h:53
Expr pow(Expr base, Expr exp)
Definition expr.cpp:235
Expr fabs(Expr x)
Definition expr.cpp:230
void print(U data, const std::string &name="", OstreamManager clout=OstreamManager(std::cout,"print"), const char delimiter=',')
bool nearZero(T a) any_platform
return true if a is close to zero
Definition util.h:402
Top level namespace for all of OpenLB.
constexpr T norm_squared(const ScalarVector< T, D, IMPL > &a) any_platform
Squared euclidean vector norm.
Vector< T, D > PhysR
Type for spatial (physical) coordinates.
std::string getName(OperatorScope scope)
Returns human-readable name of scope.
constexpr T norm(const ScalarVector< T, D, IMPL > &a) any_platform
Euclidean vector norm.
constexpr Vector< T, 3 > crossProduct3D(const ScalarVector< T, 3, IMPL > &a, const ScalarVector< T, 3, IMPL_ > &b) any_platform
Definition vector.h:263
constexpr Vector< T, D > normalize(const ScalarVector< T, D, IMPL > &a, T scale=T{1})
Definition vector.h:284
Octree - adapted from http://www.flipcode.com/archives/Octree_Implementation.shtml.
Definition of singletons: global, publicly available information.
Input in STL format – header file.
std::vector< T > getE0()
Returns Pt0-Pt1.
Definition stlReader.hh:116
std::vector< T > getE1()
Returns Pt0-Pt2.
Definition stlReader.hh:126
bool testRayIntersect(const Vector< T, 3 > &pt, const Vector< T, 3 > &dir, Vector< T, 3 > &q, T &alpha, const T &rad=T(), bool print=false)
Test intersection between ray and triangle.
Definition stlReader.hh:164
bool getPointToEdgeDistances(const Vector< T, 3 > &input, Vector< T, 3 > &output, T sensitivity=1.e-15)
Returns true if the point is on a edge (smaller than sensitivity) and gives the perpendicular distanc...
Definition stlReader.hh:328
void init()
Initializes triangle and precomputes member variables.
Definition stlReader.hh:56
bool isVortexPoint(const Vector< T, 3 > &input, Vector< T, 3 > &P, T sensitivity=1.e-15)
Returns true if is near vortex (smaller than sensitivity) and saves in P the vortex points.
Definition stlReader.hh:371
Vector< T, 3 > getCenter()
Returns center.
Definition stlReader.hh:101
bool isEdgePoint(const Vector< T, 3 > &input, Vector< T, 3 > &P1, Vector< T, 3 > &P2, T sensitivity=1.e-15)
Returns true if is near edge (smaller than sensitivity) and not near vortex and saves in P1 and P2 th...
Definition stlReader.hh:343
Vector< T, 3 > closestPtPointTriangle(const Vector< T, 3 > &pt) const
computes closest Point in a triangle to another point.
Definition stlReader.hh:271
bool isPointInside(const PhysR< T, 3 > &pt) const
Check whether a point is inside a triangle.
Definition stlReader.hh:136