00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "faceTriangulation.H"
00027 #include <OpenFOAM/plane.H>
00028
00029
00030
00031
00032 const Foam::scalar Foam::faceTriangulation::edgeRelTol = 1E-6;
00033
00034
00035
00036 Foam::label Foam::faceTriangulation::right(const label, label i)
00037 {
00038 return i;
00039 }
00040
00041
00042
00043 Foam::label Foam::faceTriangulation::left(const label size, label i)
00044 {
00045 return i ? i-1 : size-1;
00046 }
00047
00048
00049
00050
00051 Foam::tmp<Foam::vectorField> Foam::faceTriangulation::calcEdges
00052 (
00053 const face& f,
00054 const pointField& points
00055 )
00056 {
00057 tmp<vectorField> tedges(new vectorField(f.size()));
00058 vectorField& edges = tedges();
00059
00060 forAll(f, i)
00061 {
00062 point thisPt = points[f[i]];
00063 point nextPt = points[f[f.fcIndex(i)]];
00064
00065 vector vec(nextPt - thisPt);
00066 vec /= mag(vec) + VSMALL;
00067
00068 edges[i] = vec;
00069 }
00070
00071 return tedges;
00072 }
00073
00074
00075
00076 void Foam::faceTriangulation::calcHalfAngle
00077 (
00078 const vector& normal,
00079 const vector& e0,
00080 const vector& e1,
00081 scalar& cosHalfAngle,
00082 scalar& sinHalfAngle
00083 )
00084 {
00085
00086 scalar cos = max(-1, min(1, e0 & e1));
00087
00088 scalar sin = (e0 ^ e1) & normal;
00089
00090 if (sin < -ROOTVSMALL)
00091 {
00092
00093 cosHalfAngle = - Foam::sqrt(0.5*(1 + cos));
00094 sinHalfAngle = Foam::sqrt(0.5*(1 - cos));
00095 }
00096 else
00097 {
00098
00099 cosHalfAngle = Foam::sqrt(0.5*(1 + cos));
00100 sinHalfAngle = Foam::sqrt(0.5*(1 - cos));
00101 }
00102 }
00103
00104
00105
00106
00107 Foam::pointHit Foam::faceTriangulation::rayEdgeIntersect
00108 (
00109 const vector& normal,
00110 const point& rayOrigin,
00111 const vector& rayDir,
00112 const point& p1,
00113 const point& p2,
00114 scalar& posOnEdge
00115 )
00116 {
00117
00118 pointHit result(p1);
00119
00120
00121 const vector y = normal ^ rayDir;
00122
00123 posOnEdge = plane(rayOrigin, y).normalIntersect(p1, (p2-p1));
00124
00125
00126 if ((posOnEdge < 0) || (posOnEdge > 1))
00127 {
00128
00129 }
00130 else
00131 {
00132
00133 point intersectPt = p1 + posOnEdge * (p2 - p1);
00134
00135 if (((intersectPt - rayOrigin) & rayDir) < 0)
00136 {
00137
00138 }
00139 else
00140 {
00141
00142 result.setHit();
00143 result.setPoint(intersectPt);
00144 result.setDistance(mag(intersectPt - rayOrigin));
00145 }
00146 }
00147 return result;
00148 }
00149
00150
00151
00152
00153 bool Foam::faceTriangulation::triangleContainsPoint
00154 (
00155 const vector& n,
00156 const point& p0,
00157 const point& p1,
00158 const point& p2,
00159 const point& pt
00160 )
00161 {
00162 scalar area01Pt = triPointRef(p0, p1, pt).normal() & n;
00163 scalar area12Pt = triPointRef(p1, p2, pt).normal() & n;
00164 scalar area20Pt = triPointRef(p2, p0, pt).normal() & n;
00165
00166 if ((area01Pt > 0) && (area12Pt > 0) && (area20Pt > 0))
00167 {
00168 return true;
00169 }
00170 else if ((area01Pt < 0) && (area12Pt < 0) && (area20Pt < 0))
00171 {
00172 FatalErrorIn("triangleContainsPoint") << abort(FatalError);
00173 return false;
00174 }
00175 else
00176 {
00177 return false;
00178 }
00179 }
00180
00181
00182
00183
00184 void Foam::faceTriangulation::findDiagonal
00185 (
00186 const pointField& points,
00187 const face& f,
00188 const vectorField& edges,
00189 const vector& normal,
00190 const label startIndex,
00191 label& index1,
00192 label& index2
00193 )
00194 {
00195 const point& startPt = points[f[startIndex]];
00196
00197
00198 const vector& rightE = edges[right(f.size(), startIndex)];
00199 const vector leftE = -edges[left(f.size(), startIndex)];
00200
00201
00202 scalar cosHalfAngle = GREAT;
00203 scalar sinHalfAngle = GREAT;
00204 calcHalfAngle(normal, rightE, leftE, cosHalfAngle, sinHalfAngle);
00205 vector rayDir
00206 (
00207 cosHalfAngle*rightE
00208 + sinHalfAngle*(normal ^ rightE)
00209 );
00210
00211
00212 rayDir /= mag(rayDir) + VSMALL;
00213
00214
00215
00216
00217
00218
00219 label faceVertI = f.fcIndex(startIndex);
00220
00221 pointHit minInter(false, vector::zero, GREAT, true);
00222 label minIndex = -1;
00223 scalar minPosOnEdge = GREAT;
00224
00225 for (label i = 0; i < f.size() - 2; i++)
00226 {
00227 scalar posOnEdge;
00228 pointHit inter =
00229 rayEdgeIntersect
00230 (
00231 normal,
00232 startPt,
00233 rayDir,
00234 points[f[faceVertI]],
00235 points[f[f.fcIndex(faceVertI)]],
00236 posOnEdge
00237 );
00238
00239 if (inter.hit() && inter.distance() < minInter.distance())
00240 {
00241 minInter = inter;
00242 minIndex = faceVertI;
00243 minPosOnEdge = posOnEdge;
00244 }
00245
00246 faceVertI = f.fcIndex(faceVertI);
00247 }
00248
00249
00250 if (minIndex == -1)
00251 {
00252
00253
00254
00255
00256 index1 = -1;
00257 index2 = -1;
00258 return;
00259 }
00260
00261 const label leftIndex = minIndex;
00262 const label rightIndex = f.fcIndex(minIndex);
00263
00264
00265
00266
00267
00268 if (mag(minPosOnEdge) < edgeRelTol && f.fcIndex(startIndex) != leftIndex)
00269 {
00270 index1 = startIndex;
00271 index2 = leftIndex;
00272
00273 return;
00274 }
00275 if
00276 (
00277 mag(minPosOnEdge - 1) < edgeRelTol
00278 && f.fcIndex(rightIndex) != startIndex
00279 )
00280 {
00281 index1 = startIndex;
00282 index2 = rightIndex;
00283
00284 return;
00285 }
00286
00287
00288
00289
00290
00291 const point& leftPt = points[f[leftIndex]];
00292 const point& rightPt = points[f[rightIndex]];
00293
00294 minIndex = -1;
00295 scalar maxCos = -GREAT;
00296
00297
00298 faceVertI = f.fcIndex(f.fcIndex(startIndex));
00299 for (label i = 0; i < f.size() - 3; i++)
00300 {
00301 const point& pt = points[f[faceVertI]];
00302
00303 if
00304 (
00305 (faceVertI == leftIndex)
00306 || (faceVertI == rightIndex)
00307 || (triangleContainsPoint(normal, startPt, leftPt, rightPt, pt))
00308 )
00309 {
00310
00311
00312 vector edgePt0 = pt - startPt;
00313 edgePt0 /= mag(edgePt0);
00314
00315 scalar cos = rayDir & edgePt0;
00316 if (cos > maxCos)
00317 {
00318 maxCos = cos;
00319 minIndex = faceVertI;
00320 }
00321 }
00322 faceVertI = f.fcIndex(faceVertI);
00323 }
00324
00325 if (minIndex == -1)
00326 {
00327
00328
00329 index1 = startIndex;
00330
00331 if (f.fcIndex(startIndex) != leftIndex)
00332 {
00333 index2 = leftIndex;
00334 }
00335 else
00336 {
00337 index2 = rightIndex;
00338 }
00339
00340 return;
00341 }
00342
00343 index1 = startIndex;
00344 index2 = minIndex;
00345 }
00346
00347
00348
00349
00350
00351 Foam::label Foam::faceTriangulation::findStart
00352 (
00353 const face& f,
00354 const vectorField& edges,
00355 const vector& normal
00356 )
00357 {
00358 const label size = f.size();
00359
00360 scalar minCos = GREAT;
00361 label minIndex = -1;
00362
00363 forAll(f, fp)
00364 {
00365 const vector& rightEdge = edges[right(size, fp)];
00366 const vector leftEdge = -edges[left(size, fp)];
00367
00368 if (((rightEdge ^ leftEdge) & normal) < ROOTVSMALL)
00369 {
00370 scalar cos = rightEdge & leftEdge;
00371 if (cos < minCos)
00372 {
00373 minCos = cos;
00374 minIndex = fp;
00375 }
00376 }
00377 }
00378
00379 if (minIndex == -1)
00380 {
00381
00382 minCos = GREAT;
00383
00384 forAll(f, fp)
00385 {
00386 const vector& rightEdge = edges[right(size, fp)];
00387 const vector leftEdge = -edges[left(size, fp)];
00388
00389 scalar cos = rightEdge & leftEdge;
00390 if (cos < minCos)
00391 {
00392 minCos = cos;
00393 minIndex = fp;
00394 }
00395 }
00396 }
00397
00398 return minIndex;
00399 }
00400
00401
00402
00403
00404
00405
00406 bool Foam::faceTriangulation::split
00407 (
00408 const bool fallBack,
00409 const pointField& points,
00410 const face& f,
00411 const vector& normal,
00412 label& triI
00413 )
00414 {
00415 const label size = f.size();
00416
00417 if (size <= 2)
00418 {
00419 WarningIn
00420 (
00421 "split(const bool, const pointField&, const face&"
00422 ", const vector&, label&)"
00423 ) << "Illegal face:" << f
00424 << " with points " << UIndirectList<point>(points, f)()
00425 << endl;
00426
00427 return false;
00428 }
00429 else if (size == 3)
00430 {
00431
00432 triFace& tri = operator[](triI++);
00433 tri[0] = f[0];
00434 tri[1] = f[1];
00435 tri[2] = f[2];
00436
00437 return true;
00438 }
00439 else
00440 {
00441
00442
00443
00444 tmp<vectorField> tedges(calcEdges(f, points));
00445 const vectorField& edges = tedges();
00446
00447 label startIndex = findStart(f, edges, normal);
00448
00449
00450 label index1 = -1;
00451 label index2 = -1;
00452
00453 for (label iter = 0; iter < f.size(); iter++)
00454 {
00455 findDiagonal
00456 (
00457 points,
00458 f,
00459 edges,
00460 normal,
00461 startIndex,
00462 index1,
00463 index2
00464 );
00465
00466 if (index1 != -1 && index2 != -1)
00467 {
00468
00469 break;
00470 }
00471
00472
00473 startIndex = f.fcIndex(startIndex);
00474 }
00475
00476 if (index1 == -1 || index2 == -1)
00477 {
00478 if (fallBack)
00479 {
00480
00481
00482 label maxIndex = -1;
00483 scalar maxCos = -GREAT;
00484
00485 forAll(f, fp)
00486 {
00487 const vector& rightEdge = edges[right(size, fp)];
00488 const vector leftEdge = -edges[left(size, fp)];
00489
00490 scalar cos = rightEdge & leftEdge;
00491 if (cos > maxCos)
00492 {
00493 maxCos = cos;
00494 maxIndex = fp;
00495 }
00496 }
00497
00498 WarningIn
00499 (
00500 "split(const bool, const pointField&, const face&"
00501 ", const vector&, label&)"
00502 ) << "Cannot find valid diagonal on face " << f
00503 << " with points " << UIndirectList<point>(points, f)()
00504 << nl
00505 << "Returning naive triangulation starting from "
00506 << f[maxIndex] << " which might not be correct for a"
00507 << " concave or warped face" << endl;
00508
00509
00510 label fp = f.fcIndex(maxIndex);
00511
00512 for (label i = 0; i < size-2; i++)
00513 {
00514 label nextFp = f.fcIndex(fp);
00515
00516 triFace& tri = operator[](triI++);
00517 tri[0] = f[maxIndex];
00518 tri[1] = f[fp];
00519 tri[2] = f[nextFp];
00520
00521 fp = nextFp;
00522 }
00523
00524 return true;
00525 }
00526 else
00527 {
00528 WarningIn
00529 (
00530 "split(const bool, const pointField&, const face&"
00531 ", const vector&, label&)"
00532 ) << "Cannot find valid diagonal on face " << f
00533 << " with points " << UIndirectList<point>(points, f)()
00534 << nl
00535 << "Returning empty triFaceList" << endl;
00536
00537 return false;
00538 }
00539 }
00540
00541
00542
00543
00544
00545
00546
00547 label diff = 0;
00548 if (index2 > index1)
00549 {
00550 diff = index2 - index1;
00551 }
00552 else
00553 {
00554
00555 diff = index2 + size - index1;
00556 }
00557
00558 label nPoints1 = diff + 1;
00559 label nPoints2 = size - diff + 1;
00560
00561 if (nPoints1 == size || nPoints2 == size)
00562 {
00563 FatalErrorIn
00564 (
00565 "split(const bool, const pointField&, const face&"
00566 ", const vector&, label&)"
00567 ) << "Illegal split of face:" << f
00568 << " with points " << UIndirectList<point>(points, f)()
00569 << " at indices " << index1 << " and " << index2
00570 << abort(FatalError);
00571 }
00572
00573
00574
00575 face face1(nPoints1);
00576
00577 label faceVertI = index1;
00578 for (int i = 0; i < nPoints1; i++)
00579 {
00580 face1[i] = f[faceVertI];
00581 faceVertI = f.fcIndex(faceVertI);
00582 }
00583
00584
00585 face face2(nPoints2);
00586
00587 faceVertI = index2;
00588 for (int i = 0; i < nPoints2; i++)
00589 {
00590 face2[i] = f[faceVertI];
00591 faceVertI = f.fcIndex(faceVertI);
00592 }
00593
00594
00595
00596
00597
00598
00599
00600 bool splitOk =
00601 split(fallBack, points, face1, normal, triI)
00602 && split(fallBack, points, face2, normal, triI);
00603
00604
00605
00606 return splitOk;
00607 }
00608 }
00609
00610
00611
00612
00613
00614 Foam::faceTriangulation::faceTriangulation()
00615 :
00616 triFaceList()
00617 {}
00618
00619
00620
00621 Foam::faceTriangulation::faceTriangulation
00622 (
00623 const pointField& points,
00624 const face& f,
00625 const bool fallBack
00626 )
00627 :
00628 triFaceList(f.size()-2)
00629 {
00630 vector avgNormal = f.normal(points);
00631 avgNormal /= mag(avgNormal) + VSMALL;
00632
00633 label triI = 0;
00634
00635 bool valid = split(fallBack, points, f, avgNormal, triI);
00636
00637 if (!valid)
00638 {
00639 setSize(0);
00640 }
00641 }
00642
00643
00644
00645 Foam::faceTriangulation::faceTriangulation
00646 (
00647 const pointField& points,
00648 const face& f,
00649 const vector& n,
00650 const bool fallBack
00651 )
00652 :
00653 triFaceList(f.size()-2)
00654 {
00655 label triI = 0;
00656
00657 bool valid = split(fallBack, points, f, n, triI);
00658
00659 if (!valid)
00660 {
00661 setSize(0);
00662 }
00663 }
00664
00665
00666
00667 Foam::faceTriangulation::faceTriangulation(Istream& is)
00668 :
00669 triFaceList(is)
00670 {}
00671
00672
00673