OpenLB 1.7
Loading...
Searching...
No Matches
stlReader.hh
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2010-2015 Thomas Henn, Mathias J. Krause, Jonathan Jeppener-Haltenhoff
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
42// All OpenLB code is contained in this namespace.
43namespace olb {
44
45template<typename T>
47{
48 Vector<T,3> A = point[0].coords;
49 Vector<T,3> B = point[1].coords;
50 Vector<T,3> C = point[2].coords;
51 Vector<T,3> b, c;
52 T bb = 0., bc = 0., cc = 0.;
53
54 for (int i = 0; i < 3; i++) {
55 b[i] = B[i] - A[i];
56 c[i] = C[i] - A[i];
57 bb += b[i] * b[i];
58 bc += b[i] * c[i];
59 cc += c[i] * c[i];
60 }
61
62 normal[0] = b[1] * c[2] - b[2] * c[1];
63 normal[1] = b[2] * c[0] - b[0] * c[2];
64 normal[2] = b[0] * c[1] - b[1] * c[0];
65
66 T norm = util::sqrt(
67 util::pow(normal[0], 2) + util::pow(normal[1], 2) + util::pow(normal[2], 2));
68 normal[0] /= norm;
69 normal[1] /= norm;
70 normal[2] /= norm;
71
72 T D = 1.0 / (cc * bb - bc * bc);
73 T bbD = bb * D;
74 T bcD = bc * D;
75 T ccD = cc * D;
76
77 kBeta = 0.;
78 kGamma = 0.;
79 d = 0.;
80
81 for (int i = 0; i < 3; i++) {
82 uBeta[i] = b[i] * ccD - c[i] * bcD;
83 uGamma[i] = c[i] * bbD - b[i] * bcD;
84 kBeta -= A[i] * uBeta[i];
85 kGamma -= A[i] * uGamma[i];
86 d += A[i] * normal[i];
87 }
88}
89
90template<typename T>
92{
93 Vector<T,3> center( T(0) );
94
95 center[0] = (point[0].coords[0] + point[1].coords[0]
96 + point[2].coords[0]) / 3.;
97 center[1] = (point[0].coords[1] + point[1].coords[1]
98 + point[2].coords[1]) / 3.;
99 center[2] = (point[0].coords[2] + point[1].coords[2]
100 + point[2].coords[2]) / 3.;
101
102 return center;
103}
104
105template<typename T>
106std::vector<T> STLtriangle<T>::getE0()
107{
108 Vector<T,3> vec;
109 vec[0] = point[0].coords[0] - point[1].coords[0];
110 vec[1] = point[0].coords[1] - point[1].coords[1];
111 vec[2] = point[0].coords[2] - point[1].coords[2];
112 return vec;
113}
114
115template<typename T>
116std::vector<T> STLtriangle<T>::getE1()
117{
118 Vector<T,3> vec;
119 vec[0] = point[0].coords[0] - point[2].coords[0];
120 vec[1] = point[0].coords[1] - point[2].coords[1];
121 vec[2] = point[0].coords[2] - point[2].coords[2];
122 return vec;
123}
124
125template<typename T>
127{
128 // tests with T=double and T=float show that the epsilon must be increased
129 const T epsilon = std::numeric_limits<BaseType<T>>::epsilon()*T(10);
130
131 const T beta = pt * uBeta + kBeta;
132 const T gamma = pt * uGamma + kGamma;
133
134 // check if approximately equal
135 if ( util::nearZero(norm(pt - (point[0].coords + beta*(point[1].coords-point[0].coords) + gamma*(point[2].coords-point[0].coords))), epsilon) ) {
136 const T alpha = T(1) - beta - gamma;
137 return (beta >= T(0) || util::nearZero(beta, epsilon))
138 && (gamma >= T(0) || util::nearZero(gamma, epsilon))
139 && (alpha >= T(0) || util::nearZero(alpha, epsilon));
140 }
141 return false;
142}
143
144/* Schnitttest nach
145 * http://www.uninformativ.de/bin/RaytracingSchnitttests-76a577a-CC-BY.pdf
146 *
147 * Creative Commons Namensnennung 3.0 Deutschland
148 * http://creativecommons.org/licenses/by/3.0/de/
149 *
150 * P. Hofmann, 22. August 2010
151 *
152 */
153template<typename T>
155 const Vector<T,3>& dir,
156 Vector<T,3>& q, T& alpha, const T& rad,
157 bool print)
158{
159 T rn = 0.;
160 Vector<T,3> testPt = pt + rad * normal;
161 /* Vector<T,3> help; */
162
163 for (int i = 0; i < 3; i++) {
164 rn += dir[i] * normal[i];
165 }
166#ifdef OLB_DEBUG
167
168 if (print) {
169 std::cout << "Pt: " << pt[0] << " " << pt[1] << " " << pt[2] << std::endl;
170 }
171 if (print)
172 std::cout << "testPt: " << testPt[0] << " " << testPt[1] << " " << testPt[2]
173 << std::endl;
174 if (print)
175 std::cout << "PosNeg: "
176 << normal[0] * testPt[0] + normal[1] * testPt[1] + normal[2] * testPt[2]
177 - d << std::endl;
178 if (print)
179 std::cout << "Normal: " << normal[0] << " " << normal[1] << " " << normal[2]
180 << std::endl;
181#endif
182
183 // Schnitttest Flugrichtung -> Ebene
184 if (util::fabs(rn) < std::numeric_limits<T>::epsilon()) {
185#ifdef OLB_DEBUG
186 if (print) {
187 std::cout << "FALSE 1" << std::endl;
188 }
189#endif
190 return false;
191 }
192 alpha = d - testPt[0] * normal[0] - testPt[1] * normal[1] - testPt[2] * normal[2];
193 // alpha -= testPt[i] * normal[i];
194 alpha /= rn;
195
196 // Abstand Partikel Ebene
197 if (alpha < -std::numeric_limits<T>::epsilon()) {
198#ifdef OLB_DEBUG
199 if (print) {
200 std::cout << "FALSE 2" << std::endl;
201 }
202#endif
203 return false;
204 }
205 for (int i = 0; i < 3; i++) {
206 q[i] = testPt[i] + alpha * dir[i];
207 }
208 T beta = kBeta;
209 for (int i = 0; i < 3; i++) {
210 beta += uBeta[i] * q[i];
211 }
212#ifdef OLB_DEBUG
213 T dist = util::sqrt(
214 util::pow(q[0] - testPt[0], 2) + util::pow(q[1] - testPt[1], 2)
215 + util::pow(q[2] - testPt[2], 2));
216#endif
217
218 // Schnittpunkt q in der Ebene?
219 if (beta < -std::numeric_limits<T>::epsilon()) {
220#ifdef OLB_DEBUG
221
222 if (print) {
223 std::cout << "FALSE 3 BETA " << beta << " DIST " << dist << std::endl;
224 }
225#endif
226 return false;
227 }
228 T gamma = kGamma;
229 for (int i = 0; i < 3; i++) {
230 gamma += uGamma[i] * q[i];
231 }
232 if (gamma < -std::numeric_limits<T>::epsilon()) {
233#ifdef OLB_DEBUG
234 if (print) {
235 std::cout << "FALSE 4 GAMMA " << gamma << " DIST " << dist << std::endl;
236 }
237#endif
238 return false;
239 }
240 if (1. - beta - gamma < -std::numeric_limits<T>::epsilon()) {
241#ifdef OLB_DEBUG
242 if (print)
243 std::cout << "FALSE 5 VAL " << 1 - beta - gamma << " DIST " << dist
244 << std::endl;
245#endif
246 return false;
247 }
248#ifdef OLB_DEBUG
249 if (print) {
250 std::cout << "TRUE" << " GAMMA " << gamma << " BETA " << beta << std::endl;
251 }
252#endif
253 return true;
254}
255
260template<typename T>
262 const Vector<T,3>& pt) const
263{
264
265 const T nEps = -std::numeric_limits<T>::epsilon();
266 const T Eps = std::numeric_limits<T>::epsilon();
267
268 Vector<T,3> ab = point[1].coords - point[0].coords;
269 Vector<T,3> ac = point[2].coords - point[0].coords;
270 Vector<T,3> bc = point[2].coords - point[1].coords;
271
272 T snom = (pt - point[0].coords)*ab;
273 T sdenom = (pt - point[1].coords)*(point[0].coords - point[1].coords);
274
275 T tnom = (pt - point[0].coords)*ac;
276 T tdenom = (pt - point[2].coords)*(point[0].coords - point[2].coords);
277
278 if (snom < nEps && tnom < nEps) {
279 return point[0].coords;
280 }
281
282 T unom = (pt - point[1].coords)*bc;
283 T udenom = (pt - point[2].coords)*(point[1].coords - point[2].coords);
284
285 if (sdenom < nEps && unom < nEps) {
286 return point[1].coords;
287 }
288 if (tdenom < nEps && udenom < nEps) {
289 return point[2].coords;
290 }
291
292 T vc = normal*crossProduct3D(point[0].coords - pt, point[1].coords - pt);
293
294 if (vc < nEps && snom > Eps && sdenom > Eps) {
295 return point[0].coords + snom / (snom + sdenom) * ab;
296 }
297
298 T va = normal*crossProduct3D(point[1].coords - pt, point[2].coords - pt);
299
300 if (va < nEps && unom > Eps && udenom > Eps) {
301 return point[1].coords + unom / (unom + udenom) * bc;
302 }
303
304 T vb = normal*crossProduct3D(point[2].coords - pt, point[0].coords - pt);
305
306 if (vb < nEps && tnom > Eps && tdenom > Eps) {
307 return point[0].coords + tnom / (tnom + tdenom) * ac;
308 }
309
310 T u = va / (va + vb + vc);
311 T v = vb / (va + vb + vc);
312 T w = 1. - u - v;
313
314 return u * point[0].coords + v * point[1].coords + w * point[2].coords;
315}
316
317template<typename T>
318STLmesh<T>::STLmesh(std::string fName, T stlSize)
319 : _fName(fName),
320 _min(T()),
321 _max(T()),
322 _maxDist2(0),
323 clout(std::cout, "STLmesh")
324{
325 std::ifstream f(fName.c_str(), std::ios::in);
326 _triangles.reserve(10000);
327 if (!f.good()) {
328 throw std::runtime_error("STL File not valid.");
329 }
330 char buf[6];
331 buf[5] = 0;
332 f.read(buf, 5);
333 const std::string asciiHeader = "solid";
334 if (std::string(buf) == asciiHeader) {
335 f.seekg(0, std::ios::beg);
336 if (f.good()) {
337 std::string s0, s1;
338 int i = 0;
339 while (!f.eof()) {
340 f >> s0;
341 if (s0 == "facet") {
342 STLtriangle<T> tri;
343 f >> s1 >> tri.normal[0] >> tri.normal[1] >> tri.normal[2];
344 f >> s0 >> s1;
345 f >> s0 >> tri.point[0].coords[0] >> tri.point[0].coords[1]
346 >> tri.point[0].coords[2];
347 f >> s0 >> tri.point[1].coords[0] >> tri.point[1].coords[1]
348 >> tri.point[1].coords[2];
349 f >> s0 >> tri.point[2].coords[0] >> tri.point[2].coords[1]
350 >> tri.point[2].coords[2];
351 f >> s0;
352 f >> s0;
353 for (int k = 0; k < 3; k++) {
354 tri.point[0].coords[k] *= stlSize;
355 tri.point[1].coords[k] *= stlSize;
356 tri.point[2].coords[k] *= stlSize;
357 }
358 if (i == 0) {
359 _min = T();
360 _max = T();
361
362 _min[0] = tri.point[0].coords[0];
363 _min[1] = tri.point[0].coords[1];
364 _min[2] = tri.point[0].coords[2];
365
366 _max[0] = tri.point[0].coords[0];
367 _max[1] = tri.point[0].coords[1];
368 _max[2] = tri.point[0].coords[2];
369
370 _min[0] = util::min(_min[0], tri.point[1].coords[0]);
371 _min[1] = util::min(_min[1], tri.point[1].coords[1]);
372 _min[2] = util::min(_min[2], tri.point[1].coords[2]);
373
374 _max[0] = util::max(_max[0], tri.point[1].coords[0]);
375 _max[1] = util::max(_max[1], tri.point[1].coords[1]);
376 _max[2] = util::max(_max[2], tri.point[1].coords[2]);
377
378 _min[0] = util::min(_min[0], tri.point[2].coords[0]);
379 _min[1] = util::min(_min[1], tri.point[2].coords[1]);
380 _min[2] = util::min(_min[2], tri.point[2].coords[2]);
381
382 _max[0] = util::max(_max[0], tri.point[2].coords[0]);
383 _max[1] = util::max(_max[1], tri.point[2].coords[1]);
384 _max[2] = util::max(_max[2], tri.point[2].coords[2]);
385
386 }
387 else {
388 _min[0] = util::min(_min[0], tri.point[0].coords[0]);
389 _min[1] = util::min(_min[1], tri.point[0].coords[1]);
390 _min[2] = util::min(_min[2], tri.point[0].coords[2]);
391
392 _max[0] = util::max(_max[0], tri.point[0].coords[0]);
393 _max[1] = util::max(_max[1], tri.point[0].coords[1]);
394 _max[2] = util::max(_max[2], tri.point[0].coords[2]);
395
396 _min[0] = util::min(_min[0], tri.point[1].coords[0]);
397 _min[1] = util::min(_min[1], tri.point[1].coords[1]);
398 _min[2] = util::min(_min[2], tri.point[1].coords[2]);
399
400 _max[0] = util::max(_max[0], tri.point[1].coords[0]);
401 _max[1] = util::max(_max[1], tri.point[1].coords[1]);
402 _max[2] = util::max(_max[2], tri.point[1].coords[2]);
403
404 _min[0] = util::min(_min[0], tri.point[2].coords[0]);
405 _min[1] = util::min(_min[1], tri.point[2].coords[1]);
406 _min[2] = util::min(_min[2], tri.point[2].coords[2]);
407
408 _max[0] = util::max(_max[0], tri.point[2].coords[0]);
409 _max[1] = util::max(_max[1], tri.point[2].coords[1]);
410 _max[2] = util::max(_max[2], tri.point[2].coords[2]);
411 }
412
413 i++;
414 tri.init();
415 _triangles.push_back(tri);
416
417 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[1]),
418 _maxDist2);
419 _maxDist2 = util::max(distPoints(tri.point[2], tri.point[1]),
420 _maxDist2);
421 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[2]),
422 _maxDist2);
423 }
424 else if (s0 == "endsolid") {
425 break;
426 }
427 }
428 }
429 }
430 else {
431 f.close();
432 f.open(fName.c_str(), std::ios::in | std::ios::binary);
433 char comment[80];
434 f.read(comment, 80);
435
436 if (!f.good()) {
437 throw std::runtime_error("STL File not valid.");
438 }
439
440 comment[79] = 0;
441 int32_t nFacets;
442 f.read(reinterpret_cast<char *>(&nFacets), sizeof(int32_t));
443
444 if (!f.good()) {
445 throw std::runtime_error("STL File not valid.");
446 }
447
448 float v[12];
449 std::uint16_t uint16;
450 for (int32_t i = 0; i < nFacets; ++i) {
451 for (unsigned int j = 0; j < 12; ++j) {
452 f.read(reinterpret_cast<char *>(&v[j]), sizeof(float));
453 }
454 f.read(reinterpret_cast<char *>(&uint16), sizeof(std::uint16_t));
455 STLtriangle<T> tri;
456 tri.normal[0] = v[0];
457 tri.normal[1] = v[1];
458 tri.normal[2] = v[2];
459 tri.point[0].coords[0] = v[3];
460 tri.point[0].coords[1] = v[4];
461 tri.point[0].coords[2] = v[5];
462 tri.point[1].coords[0] = v[6];
463 tri.point[1].coords[1] = v[7];
464 tri.point[1].coords[2] = v[8];
465 tri.point[2].coords[0] = v[9];
466 tri.point[2].coords[1] = v[10];
467 tri.point[2].coords[2] = v[11];
468
469 for (int k = 0; k < 3; k++) {
470 tri.point[0].coords[k] *= stlSize;
471 tri.point[1].coords[k] *= stlSize;
472 tri.point[2].coords[k] *= stlSize;
473 }
474 if (i == 0) {
475 _min[0] = tri.point[0].coords[0];
476 _min[1] = tri.point[0].coords[1];
477 _min[2] = tri.point[0].coords[2];
478
479 _max[0] = tri.point[0].coords[0];
480 _max[1] = tri.point[0].coords[1];
481 _max[2] = tri.point[0].coords[2];
482
483 _min[0] = util::min(_min[0], (T) tri.point[1].coords[0]);
484 _min[1] = util::min(_min[1], (T) tri.point[1].coords[1]);
485 _min[2] = util::min(_min[2], (T) tri.point[1].coords[2]);
486
487 _max[0] = util::max(_max[0], (T) tri.point[1].coords[0]);
488 _max[1] = util::max(_max[1], (T) tri.point[1].coords[1]);
489 _max[2] = util::max(_max[2], (T) tri.point[1].coords[2]);
490
491 _min[0] = util::min(_min[0], (T) tri.point[2].coords[0]);
492 _min[1] = util::min(_min[1], (T) tri.point[2].coords[1]);
493 _min[2] = util::min(_min[2], (T) tri.point[2].coords[2]);
494
495 _max[0] = util::max(_max[0], (T) tri.point[2].coords[0]);
496 _max[1] = util::max(_max[1], (T) tri.point[2].coords[1]);
497 _max[2] = util::max(_max[2], (T) tri.point[2].coords[2]);
498
499 }
500 else {
501 _min[0] = util::min(_min[0], (T) tri.point[0].coords[0]);
502 _min[1] = util::min(_min[1], (T) tri.point[0].coords[1]);
503 _min[2] = util::min(_min[2], (T) tri.point[0].coords[2]);
504
505 _max[0] = util::max(_max[0], (T) tri.point[0].coords[0]);
506 _max[1] = util::max(_max[1], (T) tri.point[0].coords[1]);
507 _max[2] = util::max(_max[2], (T) tri.point[0].coords[2]);
508
509 _min[0] = util::min(_min[0], (T) tri.point[1].coords[0]);
510 _min[1] = util::min(_min[1], (T) tri.point[1].coords[1]);
511 _min[2] = util::min(_min[2], (T) tri.point[1].coords[2]);
512
513 _max[0] = util::max(_max[0], (T) tri.point[1].coords[0]);
514 _max[1] = util::max(_max[1], (T) tri.point[1].coords[1]);
515 _max[2] = util::max(_max[2], (T) tri.point[1].coords[2]);
516
517 _min[0] = util::min(_min[0], (T) tri.point[2].coords[0]);
518 _min[1] = util::min(_min[1], (T) tri.point[2].coords[1]);
519 _min[2] = util::min(_min[2], (T) tri.point[2].coords[2]);
520
521 _max[0] = util::max(_max[0], (T) tri.point[2].coords[0]);
522 _max[1] = util::max(_max[1], (T) tri.point[2].coords[1]);
523 _max[2] = util::max(_max[2], (T) tri.point[2].coords[2]);
524 }
525 tri.init();
526 _triangles.push_back(tri);
527
528 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[1]), _maxDist2);
529 _maxDist2 = util::max(distPoints(tri.point[2], tri.point[1]), _maxDist2);
530 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[2]), _maxDist2);
531 }
532 }
533 f.close();
534}
535
536template<typename T>
537STLmesh<T>::STLmesh(const std::vector<std::vector<T>> meshPoints, T stlSize)
538 : _fName("meshPoints.stl"),
539 _min(T()),
540 _max(T()),
541 _maxDist2(0),
542 clout(std::cout, "STLmesh")
543{
544 _triangles.reserve(10000);
545 for (size_t i = 0; i < meshPoints.size() / 3; i++) {
546 STLtriangle<T> tri;
547 tri.point[0].coords[0] = meshPoints[i*3 + 0][0];
548 tri.point[0].coords[1] = meshPoints[i*3 + 0][1];
549 tri.point[0].coords[2] = meshPoints[i*3 + 0][2];
550
551 tri.point[1].coords[0] = meshPoints[i*3 + 1][0];
552 tri.point[1].coords[1] = meshPoints[i*3 + 1][1];
553 tri.point[1].coords[2] = meshPoints[i*3 + 1][2];
554
555 tri.point[2].coords[0] = meshPoints[i*3 + 2][0];
556 tri.point[2].coords[1] = meshPoints[i*3 + 2][1];
557 tri.point[2].coords[2] = meshPoints[i*3 + 2][2];
558 for (int k = 0; k < 3; k++) {
559 tri.point[0].coords[k] *= stlSize;
560 tri.point[1].coords[k] *= stlSize;
561 tri.point[2].coords[k] *= stlSize;
562 }
563 if (i == 0) {
564 _min*=T();
565 _max*=T();
566
567 _min[0] = tri.point[0].coords[0];
568 _min[1] = tri.point[0].coords[1];
569 _min[2] = tri.point[0].coords[2];
570
571 _max[0] = tri.point[0].coords[0];
572 _max[1] = tri.point[0].coords[1];
573 _max[2] = tri.point[0].coords[2];
574
575 _min[0] = util::min(_min[0], (T) tri.point[1].coords[0]);
576 _min[1] = util::min(_min[1], (T) tri.point[1].coords[1]);
577 _min[2] = util::min(_min[2], (T) tri.point[1].coords[2]);
578
579 _max[0] = util::max(_max[0], (T) tri.point[1].coords[0]);
580 _max[1] = util::max(_max[1], (T) tri.point[1].coords[1]);
581 _max[2] = util::max(_max[2], (T) tri.point[1].coords[2]);
582
583 _min[0] = util::min(_min[0], (T) tri.point[2].coords[0]);
584 _min[1] = util::min(_min[1], (T) tri.point[2].coords[1]);
585 _min[2] = util::min(_min[2], (T) tri.point[2].coords[2]);
586
587 _max[0] = util::max(_max[0], (T) tri.point[2].coords[0]);
588 _max[1] = util::max(_max[1], (T) tri.point[2].coords[1]);
589 _max[2] = util::max(_max[2], (T) tri.point[2].coords[2]);
590
591 }
592 else {
593 _min[0] = util::min(_min[0], (T) tri.point[0].coords[0]);
594 _min[1] = util::min(_min[1], (T) tri.point[0].coords[1]);
595 _min[2] = util::min(_min[2], (T) tri.point[0].coords[2]);
596
597 _max[0] = util::max(_max[0], (T) tri.point[0].coords[0]);
598 _max[1] = util::max(_max[1], (T) tri.point[0].coords[1]);
599 _max[2] = util::max(_max[2], (T) tri.point[0].coords[2]);
600
601 _min[0] = util::min(_min[0], (T) tri.point[1].coords[0]);
602 _min[1] = util::min(_min[1], (T) tri.point[1].coords[1]);
603 _min[2] = util::min(_min[2], (T) tri.point[1].coords[2]);
604
605 _max[0] = util::max(_max[0], (T) tri.point[1].coords[0]);
606 _max[1] = util::max(_max[1], (T) tri.point[1].coords[1]);
607 _max[2] = util::max(_max[2], (T) tri.point[1].coords[2]);
608
609 _min[0] = util::min(_min[0], (T) tri.point[2].coords[0]);
610 _min[1] = util::min(_min[1], (T) tri.point[2].coords[1]);
611 _min[2] = util::min(_min[2], (T) tri.point[2].coords[2]);
612
613 _max[0] = util::max(_max[0], (T) tri.point[2].coords[0]);
614 _max[1] = util::max(_max[1], (T) tri.point[2].coords[1]);
615 _max[2] = util::max(_max[2], (T) tri.point[2].coords[2]);
616 }
617
618 tri.init();
619 _triangles.push_back(tri);
620
621 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[1]),
622 _maxDist2);
623 _maxDist2 = util::max(distPoints(tri.point[2], tri.point[1]),
624 _maxDist2);
625 _maxDist2 = util::max(distPoints(tri.point[0], tri.point[2]),
626 _maxDist2);
627 }
628}
629
630template<typename T>
632{
633 return util::pow(T(p1.coords[0] - p2.coords[0]), 2)
634 + util::pow(T(p1.coords[1] - p2.coords[1]), 2)
635 + util::pow(T(p1.coords[2] - p2.coords[2]), 2);
636}
637
638template<typename T>
639void STLmesh<T>::print(bool full)
640{
641 if (full) {
642 int i = 0;
643 clout << "Triangles: " << std::endl;
644 typename std::vector<STLtriangle<T> >::iterator it = _triangles.begin();
645
646 for (; it != _triangles.end(); ++it) {
647 clout << i++ << ": " << it->point[0].coords[0] << " " << it->point[0].coords[1]
648 << " " << it->point[0].coords[2] << " | " << it->point[1].coords[0] << " "
649 << it->point[1].coords[1] << " " << it->point[1].coords[2] << " | "
650 << it->point[2].coords[0] << " " << it->point[2].coords[1] << " "
651 << it->point[2].coords[2] << std::endl;
652 }
653 }
654 clout << "nTriangles=" << _triangles.size() << "; maxDist2=" << _maxDist2
655 << std::endl;
656 clout << "minPhysR(StlMesh)=(" << getMin()[0] << "," << getMin()[1] << ","
657 << getMin()[2] << ")";
658 clout << "; maxPhysR(StlMesh)=(" << getMax()[0] << "," << getMax()[1] << ","
659 << getMax()[2] << ")" << std::endl;
660}
661
662template<typename T>
663void STLmesh<T>::write(std::string fName)
664{
665 int rank = 0;
666#ifdef PARALLEL_MODE_MPI
667 rank = singleton::mpi().getRank();
668#endif
669 if (rank == 0) {
670 std::string fullName = singleton::directories().getVtkOutDir() + fName
671 + ".stl";
672 std::ofstream f(fullName.c_str());
673 f << "solid ascii " << fullName << "\n";
674
675 for (unsigned int i = 0; i < _triangles.size(); i++) {
676 f << "facet normal " << _triangles[i].normal[0] << " "
677 << _triangles[i].normal[1] << " " << _triangles[i].normal[2] << "\n";
678 f << " outer loop\n";
679 f << " vertex " << _triangles[i].point[0].coords[0] << " "
680 << _triangles[i].point[0].coords[1] << " " << _triangles[i].point[0].coords[2]
681 << "\n";
682 f << " vertex " << _triangles[i].point[1].coords[0] << " "
683 << _triangles[i].point[1].coords[1] << " " << _triangles[i].point[1].coords[2]
684 << "\n";
685 f << " vertex " << _triangles[i].point[2].coords[0] << " "
686 << _triangles[i].point[2].coords[1] << " " << _triangles[i].point[2].coords[2]
687 << "\n";
688 f << " endloop\n";
689 f << "endfacet\n";
690 }
691 f << "endsolid\n";
692 f.close();
693 }
694 /*if (_verbose)*/clout << "Write ... OK" << std::endl;
695}
696
697template<typename T>
698bool 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)
699{
700 std::set<unsigned int>::iterator it = tris.begin();
701 for (; it != tris.end(); ++it) {
702 if (_triangles[*it].testRayIntersect(pt, dir, q, alpha) && alpha < 1) {
703 return true;
704 }
705 }
706 return false;
707}
708
709/*
710 * STLReader functions
711 */
712template<typename T>
713STLreader<T>::STLreader(const std::string fName, T voxelSize, T stlSize,
714 int method, bool verbose, T overlap, T max)
715 : _voxelSize(voxelSize),
716 _stlSize(stlSize),
717 _overlap(overlap),
718 _fName(fName),
719 _mesh(fName, stlSize),
720 _verbose(verbose),
721 clout(std::cout, "STLreader")
722{
723 this->getName() = "STLreader";
724
725 if (_verbose) {
726 clout << "Voxelizing ..." << std::endl;
727 }
728
729 Vector<T,3> extension = _mesh.getMax() - _mesh.getMin();
730 if ( util::nearZero(max) ) {
731 max = util::max(extension[0], util::max(extension[1], extension[2])) + _voxelSize;
732 }
733 int j = 0;
734 for (; _voxelSize * util::pow(2, j) < max; j++)
735 ;
736 Vector<T,3> center;
737 T radius = _voxelSize * util::pow(2, j - 1);
738
740 for (unsigned i = 0; i < 3; i++) {
741 center[i] = (_mesh.getMin()[i] + _mesh.getMax()[i]) / 2. - _voxelSize / 4.;
742 }
743
745 _tree = new Octree<T>(center, radius, &_mesh, j, _overlap);
746
748 for (int i = 0; i < 3; i++) {
749 this->_myMin[i] = center[i] + _voxelSize / 2.;
750 this->_myMax[i] = center[i] - _voxelSize / 2.;
751 }
752 for (int i = 0; i < 3; i++) {
753 while (this->_myMin[i] > _mesh.getMin()[i]) {
754 this->_myMin[i] -= _voxelSize;
755 }
756 while (this->_myMax[i] < _mesh.getMax()[i]) {
757 this->_myMax[i] += _voxelSize;
758 }
759 this->_myMax[i] -= _voxelSize;
760 this->_myMin[i] += _voxelSize;
761 }
762
764 switch (method) {
765 case 1:
766 indicate1();
767 break;
768 case 3:
769 indicate3();
770 break;
771 case 4:
772 indicate2_Xray();
773 break;
774 case 5:
775 indicate2_Yray();
776 break;
777 default:
778 indicate2();
779 break;
780 }
781
782 if (_verbose) {
783 print();
784 }
785 if (_verbose) {
786 clout << "Voxelizing ... OK" << std::endl;
787 }
788}
789
790/*
791 * STLReader functions
792 */
793template<typename T>
794STLreader<T>::STLreader(const std::vector<std::vector<T>> meshPoints, T voxelSize, T stlSize,
795 int method, bool verbose, T overlap, T max)
796 : _voxelSize(voxelSize),
797 _stlSize(stlSize),
798 _overlap(overlap),
799 _fName("meshPoints.stl"),
800 _mesh(meshPoints, stlSize),
801 _verbose(verbose),
802 clout(std::cout, "STLreader")
803{
804 this->getName() = "STLreader";
805
806 if (_verbose) {
807 clout << "Voxelizing ..." << std::endl;
808 }
809
810 Vector<T,3> extension = _mesh.getMax() - _mesh.getMin();
811 if ( util::nearZero(max) ) {
812 max = util::max(extension[0], util::max(extension[1], extension[2])) + _voxelSize;
813 }
814 int j = 0;
815 for (; _voxelSize * util::pow(2, j) < max; j++)
816 ;
817 Vector<T,3> center;
818 T radius = _voxelSize * util::pow(2, j - 1);
819
821 for (unsigned i = 0; i < 3; i++) {
822 center[i] = (_mesh.getMin()[i] + _mesh.getMax()[i]) / 2. - _voxelSize / 4.;
823 }
824
826
827 _tree = new Octree<T>(center, radius, &_mesh, j, _overlap);
828
830 for (int i = 0; i < 3; i++) {
831 this->_myMin[i] = center[i] + _voxelSize / 2.;
832 this->_myMax[i] = center[i] - _voxelSize / 2.;
833 }
834 for (int i = 0; i < 3; i++) {
835 while (this->_myMin[i] > _mesh.getMin()[i]) {
836 this->_myMin[i] -= _voxelSize;
837 }
838 while (this->_myMax[i] < _mesh.getMax()[i]) {
839 this->_myMax[i] += _voxelSize;
840 }
841 this->_myMax[i] -= _voxelSize;
842 this->_myMin[i] += _voxelSize;
843 }
844 //automaticly choose the method with minimum extension in its direction
845
846 /*if(extension[0] == std::min_element(extension.begin(), extension.end())){
847 method = 4;
848 }
849 else if(extension[1] == std::min_element(extension.begin(), extension.end())){
850 method = 5;
851 }
852 else if(extension[2] == std::min_element(extension.begin(), extension.end())){
853 method = 0;
854 }
855 */
856
857
858 // Indicate nodes of the tree. (Inside/Outside)
859 switch (method) {
860 case 1:
861 indicate1();
862 break;
863 case 3:
864 indicate3();
865 break;
866 case 4:
867 indicate2_Xray();
868 break;
869 case 5:
870 indicate2_Yray();
871 break;
872 default:
873 indicate2();
874 break;
875 }
876
878
879 if (_verbose) {
880 print();
881 }
882 if (_verbose) {
883 clout << "Voxelizing ... OK" << std::endl;
884 }
885}
886
887template<typename T>
889{
890 delete _tree;
891}
892
893/*
894 * Old indicate function (slower, more stable)
895 * Define three rays (X-, Y-, Z-direction) for each leaf and count intersections
896 * with STL for each ray. Odd number of intersection means inside (Majority vote).
897 */
898
899template<typename T>
901{
902 std::vector<Octree<T>*> leafs;
903 _tree->getLeafs(leafs);
904 typename std::vector<Octree<T>*>::iterator it = leafs.begin();
905 Vector<T,3> dir, pt, s;
906
907 int intersections = 0;
908 int inside = 0;
909 Octree<T>* node = nullptr;
910 T step = 1. / 1000. * _voxelSize;
911 for (; it != leafs.end(); ++it) {
912 inside = 0;
913
914 pt = (*it)->getCenter();
915 intersections = 0;
916 s = pt; // + step;
917
919 dir[0] = 1;
920 dir[1] = 0;
921 dir[2] = 0;
922 while (s[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
923 node = _tree->find(s, (*it)->getMaxdepth());
924 intersections += node->testIntersection(pt, dir);
925 node->intersectRayNode(pt, dir, s);
926 s = s + step * dir;
927 }
928 inside += (intersections % 2);
929
931 intersections = 0;
932 s = pt; // + step;
933 dir[0] = 0;
934 dir[1] = 1;
935 dir[2] = 0;
936 while (s[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
937 node = _tree->find(s, (*it)->getMaxdepth());
938 intersections += node->testIntersection(pt, dir);
939 node->intersectRayNode(pt, dir, s);
940 s = s + step * dir;
941 }
942 inside += (intersections % 2);
943
945 intersections = 0;
946 s = pt; // + step;
947 dir[0] = 0;
948 dir[1] = 0;
949 dir[2] = 1;
950 while (s[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
951 node = _tree->find(s, (*it)->getMaxdepth());
952 intersections += node->testIntersection(pt, dir);
953 node->intersectRayNode(pt, dir, s);
954 s = s + step * dir;
955 }
956 inside += (intersections % 2);
957 (*it)->setInside(inside > 1);
958 }
959}
960
961/*
962 * New indicate function (faster, less stable)
963 * Define ray in Z-direction for each Voxel in XY-layer. Indicate all nodes on the fly.
964 */
965template<typename T>
966void STLreader<T>::indicate2()
967{
968 T rad = _tree->getRadius();
969 Vector<T,3> rayPt = _tree->getCenter() - rad + .5 * _voxelSize;
970 Vector<T,3> pt = rayPt;
971 Vector<T,3> rayDir;
972 rayDir[0] = 0.;
973 rayDir[1] = 0.;
974 rayDir[2] = 1.;
975 //Vector<T,3> maxEdge = _tree->getCenter() + rad;
976
977 T step = 1. / 1000. * _voxelSize;
978
979 Octree<T>* node = nullptr;
980 unsigned short rayInside = 0;
981 Vector<T,3> nodeInters;
982 while (pt[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
983 node = _tree->find(pt);
984 nodeInters = pt;
985 nodeInters[2] = node->getCenter()[2] - node->getRadius();
986 rayInside = 0;
987 while (pt[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
988 node = _tree->find(pt);
989 nodeInters = pt;
990 nodeInters[2] = node->getCenter()[2] - node->getRadius();
991 rayInside = 0;
992 while (pt[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
993 node = _tree->find(pt);
994 node->checkRay(nodeInters, rayDir, rayInside);
995 node->intersectRayNode(pt, rayDir, nodeInters);
996 pt = nodeInters + step * rayDir;
997 }
998 pt[2] = rayPt[2];
999 pt[1] += _voxelSize;
1000 }
1001 pt[1] = rayPt[1];
1002 pt[0] += _voxelSize;
1003 }
1004}
1005
1006
1007
1008/*
1009 * New indicate function (faster, less stable)
1010 * Define ray in X-direction for each Voxel in YZ-layer. Indicate all nodes on the fly.
1011 */
1012
1013template<typename T>
1014void STLreader<T>::indicate2_Xray()
1015{
1016 T rad = _tree->getRadius();
1017 Vector<T,3> rayPt = _tree->getCenter() - rad + .5 * _voxelSize;
1018 Vector<T,3> pt = rayPt;
1019 Vector<T,3> rayDir;
1020 rayDir[0] = 1.;
1021 rayDir[1] = 0.;
1022 rayDir[2] = 0.;
1023 //Vector<T,3> maxEdge = _tree->getCenter() + rad;
1024
1025 T step = 1. / 1000. * _voxelSize;
1026
1027 Octree<T>* node = nullptr;
1028 unsigned short rayInside = 0;
1029 Vector<T,3> nodeInters;
1030 while (pt[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1031 node = _tree->find(pt);
1032 nodeInters = pt;
1033 nodeInters[0] = node->getCenter()[0] - node->getRadius();
1034 rayInside = 0;
1035 while (pt[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1036 node = _tree->find(pt);
1037 nodeInters = pt;
1038 nodeInters[0] = node->getCenter()[0] - node->getRadius();
1039 rayInside = 0;
1040 while (pt[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1041 node = _tree->find(pt);
1042 node->checkRay(nodeInters, rayDir, rayInside);
1043 node->intersectRayNode(pt, rayDir, nodeInters);
1044 pt = nodeInters + step * rayDir;
1045 }
1046 pt[0] = rayPt[0];
1047 pt[1] += _voxelSize;
1048 }
1049 pt[1] = rayPt[1];
1050 pt[2] += _voxelSize;
1051 }
1052}
1053
1054/*
1055 * New indicate function (faster, less stable)
1056 * Define ray in Y-direction for each Voxel in XZ-layer. Indicate all nodes on the fly.
1057 */
1058
1059template<typename T>
1060void STLreader<T>::indicate2_Yray()
1061{
1062 T rad = _tree->getRadius();
1063 Vector<T,3> rayPt = _tree->getCenter() - rad + .5 * _voxelSize;
1064 Vector<T,3> pt = rayPt;
1065 Vector<T,3> rayDir;
1066 rayDir[0] = 0.;
1067 rayDir[1] = 1.;
1068 rayDir[2] = 0.;
1069 //Vector<T,3> maxEdge = _tree->getCenter() + rad;
1070
1071 T step = 1. / 1000. * _voxelSize;
1072
1073 Octree<T>* node = nullptr;
1074 unsigned short rayInside = 0;
1075 Vector<T,3> nodeInters;
1076 while (pt[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1077 node = _tree->find(pt);
1078 nodeInters = pt;
1079 nodeInters[1] = node->getCenter()[1] - node->getRadius();
1080 rayInside = 0;
1081 while (pt[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1082 node = _tree->find(pt);
1083 nodeInters = pt;
1084 nodeInters[1] = node->getCenter()[1] - node->getRadius();
1085 rayInside = 0;
1086 while (pt[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1087 node = _tree->find(pt);
1088 node->checkRay(nodeInters, rayDir, rayInside);
1089 node->intersectRayNode(pt, rayDir, nodeInters);
1090 pt = nodeInters + step * rayDir;
1091 }
1092 pt[1] = rayPt[1];
1093 pt[0] += _voxelSize;
1094 }
1095 pt[0] = rayPt[0];
1096 pt[2] += _voxelSize;
1097 }
1098}
1099
1100/*
1101 * Double ray approach: two times (X-, Y-, Z-direction) for each leaf.
1102 * Could be use to deal with double layer triangles and face intersections.
1103 */
1104template<typename T>
1105void STLreader<T>::indicate3()
1106{
1107 std::vector<Octree<T>*> leafs;
1108 _tree->getLeafs(leafs);
1109 typename std::vector<Octree<T>*>::iterator it = leafs.begin();
1110
1111 Vector<T,3> dir, pt, s;
1112 Octree<T>* node = nullptr;
1113 T step = 1. / 1000. * _voxelSize;
1114 int intersections;
1115 int sum_intersections;
1116
1117 for (; it != leafs.end(); ++it) {
1118 pt = (*it)->getCenter();
1119 intersections = 0;
1120 sum_intersections = 0;
1121 s = pt; // + step;
1122
1124 dir[0] = 1;
1125 dir[1] = 0;
1126 dir[2] = 0;
1127 while (s[0] < _mesh.getMax()[0] + std::numeric_limits<T>::epsilon()) {
1128 node = _tree->find(s, (*it)->getMaxdepth());
1129 intersections = node->testIntersection(pt, dir);
1130 node->intersectRayNode(pt, dir, s);
1131 s = s + step * dir;
1132 if (intersections > 0) {
1133 sum_intersections++;
1134 break;
1135 }
1136 }
1137
1139 intersections = 0;
1140 s = pt; // + step;
1141 dir[0] = 0;
1142 dir[1] = 1;
1143 dir[2] = 0;
1144 while (s[1] < _mesh.getMax()[1] + std::numeric_limits<T>::epsilon()) {
1145 node = _tree->find(s, (*it)->getMaxdepth());
1146 intersections = node->testIntersection(pt, dir);
1147 node->intersectRayNode(pt, dir, s);
1148 s = s + step * dir;
1149 if (intersections > 0) {
1150 sum_intersections++;
1151 break;
1152 }
1153 }
1154
1156 intersections = 0;
1157 s = pt; // + step;
1158 dir[0] = 0;
1159 dir[1] = 0;
1160 dir[2] = 1;
1161 while (s[2] < _mesh.getMax()[2] + std::numeric_limits<T>::epsilon()) {
1162 node = _tree->find(s, (*it)->getMaxdepth());
1163 intersections = node->testIntersection(pt, dir);
1164 node->intersectRayNode(pt, dir, s);
1165 s = s + step * dir;
1166 if (intersections > 0) {
1167 sum_intersections++;
1168 break;
1169 }
1170 }
1171
1173 intersections = 0;
1174 s = pt; // + step;
1175 dir[0] = -1;
1176 dir[1] = 0;
1177 dir[2] = 0;
1178 while (s[0] > _mesh.getMin()[0] - std::numeric_limits<T>::epsilon()) {
1179 node = _tree->find(s, (*it)->getMaxdepth());
1180 intersections = node->testIntersection(pt, dir);
1181 node->intersectRayNode(pt, dir, s);
1182 s = s + step * dir;
1183 if (intersections > 0) {
1184 sum_intersections++;
1185 break;
1186 }
1187 }
1188
1190 intersections = 0;
1191 s = pt; // + step;
1192 dir[0] = 0;
1193 dir[1] = -1;
1194 dir[2] = 0;
1195 while (s[1] > _mesh.getMin()[1] - std::numeric_limits<T>::epsilon()) {
1196 node = _tree->find(s, (*it)->getMaxdepth());
1197 intersections = node->testIntersection(pt, dir);
1198 node->intersectRayNode(pt, dir, s);
1199 s = s + step * dir;
1200 if (intersections > 0) {
1201 sum_intersections++;
1202 break;
1203 }
1204 }
1205
1207 intersections = 0;
1208 s = pt; // + step;
1209 dir[0] = 0;
1210 dir[1] = 0;
1211 dir[2] = -1;
1212 while (s[2] > _mesh.getMin()[2] - 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 (*it)->setInside(sum_intersections > 5);
1223 }
1224}
1225
1226template<typename T>
1227bool STLreader<T>::operator() (bool output[], const T input[])
1228{
1229 output[0] = false;
1230 T coords = _tree->getRadius();
1231 Vector<T,3> c(_tree->getCenter());
1232 if (c[0] - coords < input[0] && input[0] < c[0] + coords && c[1] - coords < input[1]
1233 && input[1] < c[1] + coords && c[2] - coords < input[2] && input[2] < c[2] + coords) {
1234 std::vector<T> tmp(input, input + 3);
1235 output[0] = _tree->find(tmp)->getInside();
1236 }
1237 return true;
1238}
1239
1240template<typename T>
1241bool STLreader<T>::distance(T& distance, const Vector<T,3>& origin,
1242 const Vector<T,3>& direction, int iC)
1243{
1244 Octree<T>* node = nullptr;
1245 Vector<T,3> dir(direction);
1246 dir = normalize(dir);
1247 Vector<T,3> extends = _mesh.getMax() - _mesh.getMin();
1248 Vector<T,3> pt(origin);
1249 Vector<T,3> q;
1250 Vector<T,3> s;
1251 Vector<T,3> center = _mesh.getMin() + 1 / 2. * extends;
1252 T step = _voxelSize / 1000., a = 0;
1253
1254 for (int i = 0; i < 3; i++) {
1255 extends[i] /= 2.;
1256 }
1257
1258 if (!(_mesh.getMin()[0] < origin[0] && origin[0] < _mesh.getMax()[0]
1259 && _mesh.getMin()[1] < origin[1] && origin[1] < _mesh.getMax()[1]
1260 && _mesh.getMin()[2] < origin[2] && origin[2] < _mesh.getMax()[2])) {
1261 T t = T(), d = T();
1262 bool foundQ = false;
1263
1264 if (dir[0] > 0) {
1265 d = _mesh.getMin()[0];
1266 t = (d - origin[0]) / dir[0];
1267 pt[0] = origin[0] + (t + step) * dir[0];
1268 pt[1] = origin[1] + (t + step) * dir[1];
1269 pt[2] = origin[2] + (t + step) * dir[2];
1270
1271 if (_mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]
1272 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1273 foundQ = true;
1274 }
1275 }
1276 else if (dir[0] < 0) {
1277 d = _mesh.getMax()[0];
1278 t = (d - origin[0]) / dir[0];
1279 pt[0] = origin[0] + (t + step) * dir[0];
1280 pt[1] = origin[1] + (t + step) * dir[1];
1281 pt[2] = origin[2] + (t + step) * dir[2];
1282 if (_mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]
1283 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1284 foundQ = true;
1285 }
1286 }
1287
1288 if (dir[1] > 0 && !foundQ) {
1289 d = _mesh.getMin()[1];
1290 t = (d - origin[1]) / dir[1];
1291 pt[0] = origin[0] + (t + step) * dir[0];
1292 pt[1] = origin[1] + (t + step) * dir[1];
1293 pt[2] = origin[2] + (t + step) * dir[2];
1294 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1295 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1296 foundQ = true;
1297 }
1298 }
1299 else if (dir[1] < 0 && !foundQ) {
1300 d = _mesh.getMax()[1];
1301 t = (d - origin[1]) / dir[1];
1302 pt[0] = origin[0] + (t + step) * dir[0];
1303 pt[1] = origin[1] + (t + step) * dir[1];
1304 pt[2] = origin[2] + (t + step) * dir[2];
1305 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1306 && _mesh.getMin()[2] < pt[2] && pt[2] < _mesh.getMax()[2]) {
1307 foundQ = true;
1308 }
1309 }
1310
1311 if (dir[2] > 0 && !foundQ) {
1312 d = _mesh.getMin()[2];
1313 t = (d - origin[2]) / dir[2];
1314 pt[0] = origin[0] + (t + step) * dir[0];
1315 pt[1] = origin[1] + (t + step) * dir[1];
1316 pt[2] = origin[2] + (t + step) * dir[2];
1317 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1318 && _mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]) {
1319 foundQ = true;
1320 }
1321 }
1322 else if (dir[2] < 0 && !foundQ) {
1323 d = _mesh.getMax()[2];
1324 t = (d - origin[2]) / dir[2];
1325 pt[0] = origin[0] + (t + step) * dir[0];
1326 pt[1] = origin[1] + (t + step) * dir[1];
1327 pt[2] = origin[2] + (t + step) * dir[2];
1328 if (_mesh.getMin()[0] < pt[0] && pt[0] < _mesh.getMax()[0]
1329 && _mesh.getMin()[1] < pt[1] && pt[1] < _mesh.getMax()[1]) {
1330 foundQ = true;
1331 }
1332 }
1333
1334 if (!foundQ) {
1335 return false;
1336 }
1337 }
1338
1339 while ((util::fabs(pt[0] - center[0]) < extends[0])
1340 && (util::fabs(pt[1] - center[1]) < extends[1])
1341 && (util::fabs(pt[2] - center[2]) < extends[2])) {
1342 node = _tree->find(pt);
1343 if (node->closestIntersection(Vector<T,3>(origin), dir, q, a)) {
1344 Vector<T,3> vek(q - Vector<T,3>(origin));
1345 distance = norm(vek);
1346 return true;
1347 }
1348 else {
1349 Octree<T>* tmpNode = _tree->find(pt);
1350 tmpNode->intersectRayNode(pt, dir, s);
1351 for (int i = 0; i < 3; i++) {
1352 pt[i] = s[i] + step * dir[i];
1353 }
1354 }
1355 }
1356
1357 return false;
1358}
1359
1360template<typename T>
1362{
1363 // Check if the position is on the corner of a triangle
1364 unsigned countTriangles = 0;
1365 Vector<T,3> normal(T(0));
1366
1367 // TODO: Calculate angle-weighted psuedonormal (see 10.1109/TVCG.2005.49) in case the point lies on corners
1368 // Edges correspond to an unweighted average (as calculated below) anyway
1369 for (const STLtriangle<T>& triangle : _mesh.getTriangles()) {
1370 if (triangle.isPointInside(pt)) {
1371 ++countTriangles;
1372 normal+=triangle.getNormal();
1373 }
1374 }
1375 if (countTriangles > 0) {
1376 return normal / countTriangles;
1377 }
1378
1379 // If the provided point isn't located on the surface, return the predefined fallback
1380 return fallbackNormal;
1381}
1382
1383template<typename T>
1384Vector<T,3> STLreader<T>::evalSurfaceNormal(const Vector<T,3>& origin)
1385{
1386 Vector<T,3> normal(0.);
1387 Vector<T,3> closestPointOnSurface(0.);
1388 const STLtriangle<T>* closestTriangle = nullptr;
1389 T distance = std::numeric_limits<T>::max();
1390
1391 for (const STLtriangle<T>& triangle : _mesh.getTriangles()) {
1392 PhysR<T,3> const pointOnTriangle = triangle.closestPtPointTriangle(origin);
1393 PhysR<T,3> const currDistance = origin - pointOnTriangle;
1394 T currDistanceNorm = norm_squared(currDistance);
1395 if (distance > currDistanceNorm) {
1396 normal = currDistance;
1397 distance = currDistanceNorm;
1398 closestPointOnSurface = pointOnTriangle;
1399 closestTriangle = &triangle;
1400 }
1401 }
1402
1403 distance = util::sqrt(distance);
1404 if (!util::nearZero(distance)) {
1405 const short signAtInput = evalSignForSignedDistance(origin);
1406 normal = normal/distance;
1407 normal *= signAtInput;
1408 }
1409 else {
1410 normal = evalNormalOnSurface(closestPointOnSurface, closestTriangle->getNormal());
1411 // The below isn't necessary because all normals are set to point outside by default (see constructor)
1412
1413 /* const T tmpStepSize = _voxelSize * T{0.01}; */
1414 /* const short signBefore = evalSignForSignedDistance(origin - tmpStepSize * normal, *closestTriangle); */
1415 /* const short signAfter = evalSignForSignedDistance(origin + tmpStepSize * normal, *closestTriangle); */
1416
1417 /* if(signBefore > 0 || signAfter < 0) { */
1418 /* normal *= T{-1}; */
1419 /* } */
1420 }
1421
1422 return normal;
1423}
1424
1425template<typename T>
1426short STLreader<T>::evalSignForSignedDistance(const Vector<T,3>& pseudonormal, const Vector<T,3>& distance)
1427{
1428 const T projection = pseudonormal * distance;
1429
1430 if(projection > 0) {
1431 return 1;
1432 }
1433 else if (projection < 0 ) {
1434 return -1;
1435 }
1436 return 0;
1437}
1438
1439template<typename T>
1440short STLreader<T>::evalSignForSignedDistance(const Vector<T,3>& pt)
1441{
1442 T windingNumber{0};
1443
1444 for (const STLtriangle<T>& triangle : _mesh.getTriangles()) {
1445 const PhysR<T,3> a = triangle.point[0].coords - pt;
1446 const PhysR<T,3> b = triangle.point[1].coords - pt;
1447 const PhysR<T,3> c = triangle.point[2].coords - pt;
1448
1449 const T aNorm = norm(a);
1450 const T bNorm = norm(b);
1451 const T cNorm = norm(c);
1452
1453 const T numerator = a * crossProduct3D(b, c);
1454 const T denominator = aNorm * bNorm * cNorm + (a*b) * cNorm + (b*c) * aNorm + (c*a) * bNorm;
1455
1456 windingNumber += util::atan2(numerator,denominator);
1457 }
1458 windingNumber *= 2;
1459
1460 if(int(util::round(windingNumber))%2 == 0) {
1461 return 1;
1462 }
1463 return -1;
1464}
1465
1466template<typename T>
1468{
1469 Vector<T,3> distance;
1470 Vector<T,3> closestPointOnSurface;
1471 T distanceNorm = std::numeric_limits<T>::max();
1472
1473 for (const STLtriangle<T>& triangle : _mesh.getTriangles()) {
1474 const PhysR<T,3> currPointOnTriangle = triangle.closestPtPointTriangle(input);
1475 const Vector<T,3> currDistance = input - currPointOnTriangle;
1476 T currDistanceNorm = norm_squared(currDistance);
1477 if (distanceNorm > currDistanceNorm) {
1478 distanceNorm = currDistanceNorm;
1479 distance = currDistance;
1480 closestPointOnSurface = currPointOnTriangle;
1481 }
1482 }
1483
1484 return util::sqrt(distanceNorm) * evalSignForSignedDistance(input);
1485}
1486
1487template <typename T>
1489{
1490 return evalSurfaceNormal(pos);
1491}
1492
1493template<typename T>
1495{
1496 _mesh.print();
1497 _tree->print();
1498 clout << "voxelSize=" << _voxelSize << "; stlSize=" << _stlSize << std::endl;
1499 clout << "minPhysR(VoxelMesh)=(" << this->_myMin[0] << "," << this->_myMin[1]
1500 << "," << this->_myMin[2] << ")";
1501 clout << "; maxPhysR(VoxelMesh)=(" << this->_myMax[0] << ","
1502 << this->_myMax[1] << "," << this->_myMax[2] << ")" << std::endl;
1503}
1504
1505template<typename T>
1507{
1508 _tree->write(_fName);
1509}
1510
1511template<typename T>
1512void STLreader<T>::writeSTL(std::string stlName)
1513{
1514 if (stlName == "") {
1515 _mesh.write(_fName);
1516 }
1517 else {
1518 _mesh.write(stlName);
1519 }
1520}
1521
1522template<typename T>
1524{
1525 unsigned int noTris = _mesh.triangleSize();
1526 Vector<T,3> center;
1527 //Octree<T>* node = nullptr;
1528 for (unsigned int i = 0; i < noTris; i++) {
1529 center = _mesh.getTri(i).getCenter();
1530 if (_tree->find(
1531 center + _mesh.getTri(i).normal * util::sqrt(3.) * _voxelSize)->getInside()) {
1532 // cout << "Wrong direction" << std::endl;
1533 Vector<T,3> pt(_mesh.getTri(i).point[0].coords);
1534 _mesh.getTri(i).point[0].coords = _mesh.getTri(i).point[2].coords;
1535 _mesh.getTri(i).point[2].coords = pt;
1536 _mesh.getTri(i).init();
1537 // _mesh.getTri(i).getNormal()[0] *= -1.;
1538 // _mesh.getTri(i).getNormal()[1] *= -1.;
1539 // _mesh.getTri(i).getNormal()[2] *= -1.;
1540 }
1541 }
1542}
1543
1544template<typename T>
1546{
1547 std::vector<Octree<T>*> leafs;
1548 _tree->getLeafs(leafs);
1549 for (auto it = leafs.begin(); it != leafs.end(); ++it) {
1550 if ((*it)->getBoundaryNode()) {
1551 (*it)->setInside(true);
1552 }
1553 }
1554}
1555
1556} // namespace olb
1557
1558#endif
std::string & getName()
read and write access to name
Definition genericF.hh:51
bool closestIntersection(const Vector< T, 3 > &pt, const Vector< T, 3 > &direction, Vector< T, 3 > &q, T &a)
Test intersection of ray with all triangles in Octree q contains point of closest intersection to pt ...
Definition octree.hh:669
void intersectRayNode(const Vector< T, 3 > &pt, const Vector< T, 3 > &dir, Vector< T, 3 > &s)
Computes intersection of ray with Octree boundaries.
Definition octree.hh:676
int testIntersection(const Vector< T, 3 > &pt, const Vector< T, 3 > &dir, bool print=false)
Test intersection of ray with all triangles in Octree returns number of intersections.
Definition octree.hh:324
Octree< T > * find(const Vector< T, 3 > &, const int &maxDepth=0)
Find the node containing the first param with remaining maxDepth.
Definition octree.hh:263
STLmesh(std::string, T stlSize=1.)
Constructs a new STLmesh from a file.
Definition stlReader.hh:318
bool testRayIntersect(const std::set< unsigned int > &tris, const Vector< T, 3 > &pt, const Vector< T, 3 > &dir, Vector< T, 3 > &q, T &alpha)
Compute intersection between Ray and set of triangles; returns true if intersection is found.
Definition stlReader.hh:698
void print(bool full=false)
Prints console output.
Definition stlReader.hh:639
void write(std::string fName)
Writes STL mesh in Si units.
Definition stlReader.hh:663
STLreader(const std::string fName, T voxelSize, T stlSize=1, int method=2, bool verbose=false, T overlap=0., T max=0.)
Constructs a new STLreader from a file.
Definition stlReader.hh:713
T signedDistance(const Vector< T, 3 > &input) override
Computes signed distance to closest triangle in direction of the surface normal.
void print()
Prints console output.
void writeSTL(std::string stlName="")
Writes STL mesh in Si units.
void setNormalsOutside()
Rearranges normals of triangles to point outside of geometry.
Vector< T, 3 > surfaceNormal(const Vector< T, 3 > &pos, const T meshSize=0) override
Finds and returns normal of the closest surface (triangle)
~STLreader() override
Definition stlReader.hh:888
void setBoundaryInsideNodes()
Every octree leaf intersected by the STL will be part of the inside nodes.
bool operator()(bool output[], const T input[]) override
Returns whether node is inside or not.
bool distance(T &distance, const Vector< T, 3 > &origin, const Vector< T, 3 > &direction, int iC=-1) override
Computes distance to closest triangle intersection.
void writeOctree()
Writes Octree.
Plain old scalar vector.
Definition vector.h:47
std::string getVtkOutDir() const
Definition singleton.h:97
int getRank() const
Returns the process ID.
Wrapper functions that simplify the use of MPI.
platform_constant int c[Q][D]
platform_constant Fraction s[Q]
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:143
MpiManager & mpi()
Directories & directories()
Definition singleton.h:150
cpu::simd::Pack< T > sqrt(cpu::simd::Pack< T > value)
Definition pack.h:100
bool distance(S &distance, const Vector< S, D > &origin, const Vector< S, D > &direction, S precision, S pitch, F1 isInside, F2 isInsideBoundingBox)
cpu::simd::Pack< T > min(cpu::simd::Pack< T > rhs, cpu::simd::Pack< T > lhs)
Definition pack.h:124
ADf< T, DIM > atan2(const T &y, const ADf< T, DIM > &x)
Definition aDiff.h:623
ADf< T, DIM > round(const ADf< T, DIM > &a)
Definition aDiff.h:928
cpu::simd::Pack< T > max(cpu::simd::Pack< T > rhs, cpu::simd::Pack< T > lhs)
Definition pack.h:130
cpu::simd::Pack< T > pow(cpu::simd::Pack< T > base, cpu::simd::Pack< T > exp)
Definition pack.h:112
cpu::simd::Pack< T > fabs(cpu::simd::Pack< T > value)
Definition pack.h:106
bool nearZero(const ADf< T, DIM > &a)
Definition aDiff.h:1087
Top level namespace for all of OpenLB.
constexpr T norm(const ScalarVector< T, D, IMPL > &a)
Euclidean vector norm.
constexpr T norm_squared(const ScalarVector< T, D, IMPL > &a)
Squared 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:224
constexpr Vector< T, D > normalize(const ScalarVector< T, D, IMPL > &a, T scale=T{1})
Definition vector.h:245
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:106
std::vector< T > getE1()
Returns Pt0-Pt2.
Definition stlReader.hh:116
Vector< T, 3 > normal
normal of triangle
Definition stlReader.h:95
std::vector< Vertex< T, 3 > > point
A triangle contains 3 Points.
Definition stlReader.h:92
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:154
void init()
Initializes triangle and precomputes member variables.
Definition stlReader.hh:46
Vector< T, 3 > getCenter()
Returns center.
Definition stlReader.hh:91
Vector< T, 3 > closestPtPointTriangle(const Vector< T, 3 > &pt) const
computes closest Point in a triangle to another point.
Definition stlReader.hh:261
bool isPointInside(const PhysR< T, 3 > &pt) const
Check whether a point is inside a triangle.
Definition stlReader.hh:126
Vector< T, D > coords
Point coordinates in SI units.
Definition stlReader.h:68