triSurfaceTools.H 21.2 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
OpenFOAM bot's avatar
OpenFOAM bot committed
6 7
     \\/     M anipulation  |
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8 9
    Copyright (C) 2011-2018 OpenFOAM Foundation
    Copyright (C) 2019 OpenCFD Ltd.
10 11 12 13
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

14 15 16 17
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
18 19 20 21 22 23 24

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
25
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
26 27 28 29 30

Class
    Foam::triSurfaceTools

Description
31
    A collection of tools for triSurface.
32

33 34 35 36 37 38 39 40 41 42 43
Note
    The curvature calculation is an implementation of the algorithm from:
    \verbatim
        "Estimating Curvatures and their Derivatives on Triangle Meshes"
        by S. Rusinkiewicz
        3DPVT'04 Proceedings of the 3D Data Processing,
        Visualization, and Transmission, 2nd International Symposium
        Pages 486-493
        http://gfx.cs.princeton.edu/pubs/_2004_ECA/curvpaper.pdf
    \endverbatim

44 45
SourceFiles
    triSurfaceTools.C
46 47
    triSurfaceCloseness.C
    triSurfaceCurvature.C
48 49 50 51 52 53 54 55

\*---------------------------------------------------------------------------*/

#ifndef triSurfaceTools_H
#define triSurfaceTools_H

#include "boolList.H"
#include "pointField.H"
56 57
#include "vectorField.H"
#include "triadFieldFwd.H"
58
#include "DynamicList.H"
59
#include "HashSet.H"
60
#include "Map.H"
61
#include "FixedList.H"
62
#include "Pair.H"
63 64 65 66 67 68 69 70 71 72
#include "vector2D.H"
#include "triPointRef.H"
#include "surfaceLocation.H"

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{

// Forward declaration of classes
73
class boundBox;
74 75 76 77
class edge;
class labelledTri;
class polyBoundaryMesh;
class plane;
78 79
class triSurface;
class face;
80
class Time;
81 82
template<class Face> class MeshedSurface;

83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101
/*---------------------------------------------------------------------------*\
                           Class triSurfaceTools Declaration
\*---------------------------------------------------------------------------*/

class triSurfaceTools
{
    // Private Member Functions

        // Refinement

            enum refineType
            {
                NONE,
                RED,
                GREEN
            };
            static void calcRefineStatus
            (
                const triSurface& surf,
102
                const label facei,
103 104 105 106 107
                List<refineType>& refine
            );
            static void greenRefine
            (
                const triSurface& surf,
108
                const label facei,
109
                const label edgeI,
110
                const label newPointi,
111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136
                DynamicList<labelledTri>& newFaces
            );
            static triSurface doRefine
            (
                const triSurface& surf,
                const List<refineType>& refineStatus
            );


        // Coarsening

            static scalar faceCosAngle
            (
                const point& pStart,
                const point& pEnd,
                const point& pLeft,
                const point& pRight
            );

            static void protectNeighbours
            (
                const triSurface& surf,
                const label vertI,
                labelList& faceStatus
            );

137
            //- Faces to collapse because of edge collapse
138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165
            static labelHashSet getCollapsedFaces
            (
                const triSurface& surf,
                label edgeI
            );

            // Return value of faceUsed for faces using vertI (local numbering).
            // Used internally.
            static label vertexUsesFace
            (
                const triSurface& surf,
                const labelHashSet& faceUsed,
                const label vertI
            );

            // Get new connections between faces (because of edge collapse)
            // in form of tables:
            //  - given edge get other edge
            //  - given edge get other face
            // A face using point v1 on edge will get connected to a face using
            //  point v2 if they share a common vertex
            //  (but not a common edge since then the triangles collapse to
            //  nothing)
            static void getMergedEdges
            (
                const triSurface& surf,
                const label edgeI,
                const labelHashSet& collapsedFaces,
166 167
                Map<label>& edgeToEdge,
                Map<label>& edgeToFace
168 169
            );

170
            //- Calculates (cos of) angle across edgeI of facei,
171 172 173 174 175 176 177 178
            //  taking into account updated addressing (resulting from edge
            //  collapse)
            static scalar edgeCosAngle
            (
                const triSurface& surf,
                const label v1,
                const point& pt,
                const labelHashSet& collapsedFaces,
179 180
                const Map<label>& edgeToEdge,
                const Map<label>& edgeToFace,
181
                const label facei,
182 183 184 185
                const label edgeI
            );

            //- Calculate minimum (cos of) edge angle using addressing from
Henry Weller's avatar
Henry Weller committed
186
            //  collapsing
187 188 189 190 191 192 193 194
            //  edge to v1 at pt. Returns 1 if v1 is on edge without neighbours
            //  (and hence no edge angle can be defined)
            static scalar collapseMinCosAngle
            (
                const triSurface& surf,
                const label v1,
                const point& pt,
                const labelHashSet& collapsedFaces,
195 196
                const Map<label>& edgeToEdge,
                const Map<label>& edgeToFace
197 198 199 200 201 202 203 204 205
            );

            //- Like collapseMinCosAngle but return true for value < minCos
            bool collapseCreatesFold
            (
                const triSurface& surf,
                const label v1,
                const point& pt,
                const labelHashSet& collapsedFaces,
206 207
                const Map<label>& edgeToEdge,
                const Map<label>& edgeToFace,
208 209 210 211 212 213 214 215 216
                const scalar minCos
            );

            ////- Checks if edge collapse creates triangles on top of each
            ////  other
            //static bool collapseCreatesDuplicates
            //(
            //    const triSurface& surf,
            //    const label edgeI,
217
            //    const Map<bool>& collapsedFaces
218 219 220 221 222 223 224 225 226 227 228 229 230 231
            //);

        // Tracking

            //- Finds the triangle edge/point cut by the plane between
            //  a point inside/on edge of a triangle and a point outside.
            //  Returns
            //  - location on edge/point and hit()
            //  - or miss() if no intersection found
            static surfaceLocation cutEdge
            (
                const triSurface& s,
                const label triI,
                const label excludeEdgeI,
232
                const label excludePointi,
233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260
                const point& triPoint,
                const plane& cutPlane,
                const point& toPoint
            );

            //- Checks if current is on the same triangle as the endpoint
            //  and shifts it there. If so updates current and sets a hit.
            static void snapToEnd
            (
                const triSurface& s,
                const surfaceLocation& endInfo,
                surfaceLocation& current
            );

            //- Visits faces eFaces around start. Does not visit triangle
            //  start.triangle() nor edge excludeEdgeI.
            //  Returns edge, triangle (if more than one choice) which gets
            //  us nearer endpoint.
            //  Returns
            //  - hit() if triangle contains endpoint
            //  - triangle()=-1 if no triangle found
            //  - nearest triangle/edge otherwise
            static surfaceLocation visitFaces
            (
                const triSurface& s,
                const labelList& eFaces,
                const surfaceLocation& start,
                const label excludeEdgeI,
261
                const label excludePointi,
262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303
                const surfaceLocation& end,
                const plane& cutPlane
            );


public:

    // OBJ writing

        //- Write pointField to OBJ format file
        static void writeOBJ
        (
            const fileName& fName,
            const pointField& pts
        );

        //- Write vertex subset to OBJ format file
        static void writeOBJ
        (
            const triSurface& surf,
            const fileName& fName,
            const boolList& markedVerts
        );


    // Additional addressing

        //- Get all triangles using edge endpoint
        static void getVertexTriangles
        (
            const triSurface& surf,
            const label edgeI,
            labelList& edgeTris
        );

        //- Get all vertices (local numbering) connected to vertices of edge
        static labelList getVertexVertices
        (
            const triSurface& surf,
            const edge& e
        );

304
        //- Get face connected to edge not facei
305 306 307
        static label otherFace
        (
            const triSurface& surf,
308
            const label facei,
309 310 311
            const label edgeI
        );

312
        //- Get the two edges on facei counterclockwise after edgeI
313 314 315
        static void otherEdges
        (
            const triSurface& surf,
316
            const label facei,
317 318 319 320 321
            const label edgeI,
            label& e1,
            label& e2
        );

322
        //- Get the two vertices (local numbering) on facei counterclockwise
323 324 325 326
        //  vertI
        static void otherVertices
        (
            const triSurface& surf,
327
            const label facei,
328 329 330 331 332 333 334 335 336
            const label vertI,
            label& vert1I,
            label& vert2I
        );

        //- Get edge opposite vertex (local numbering)
        static label oppositeEdge
        (
            const triSurface& surf,
337
            const label facei,
338 339 340 341 342 343 344
            const label vertI
        );

        //- Get vertex (local numbering) opposite edge
        static label oppositeVertex
        (
            const triSurface& surf,
345
            const label facei,
346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445
            const label edgeI
        );

        //- Returns edge label connecting v1, v2 (local numbering)
        static label getEdge
        (
            const triSurface& surf,
            const label vert1I,
            const label vert2I
        );

        //- Return index of triangle (or -1) using all three edges
        static label getTriangle
        (
            const triSurface& surf,
            const label e0I,
            const label e1I,
            const label e2I
        );

    // Coarsening

        //- Create new triSurface by collapsing edges to edge mids.
        static triSurface collapseEdges
        (
            const triSurface& surf,
            const labelList& collapsableEdges
        );


        //- Face collapse status.
        //  anyEdge: any edge can be collapsed
        //  noEdge: no edge can be collapsed
        //  collapsed: already collapsed
        //  >0: edge label that can be collapsed
        static const label ANYEDGE;
        static const label NOEDGE;
        static const label COLLAPSED;

        //- Create new triSurface by collapsing edges to specified
        //  positions. faceStatus allows
        //  explicit control over which faces need to be protected (see above).
        //  faceStatus gets updated to protect collapsing already collapsed
        //  faces.
        static triSurface collapseEdges
        (
            const triSurface& surf,
            const labelList& collapsableEdges,
            const pointField& edgeMids,
            labelList& faceStatus
        );


    // Refinement

        //- Refine edges by splitting to opposite vertex
        static triSurface greenRefine
        (
            const triSurface& surf,
            const labelList& refineEdges
        );

        //- Refine face by splitting all edges. Neighbouring face is
        //  greenRefine'd.
        static triSurface redGreenRefine
        (
            const triSurface& surf,
            const labelList& refineFaces
        );


    // Geometric

        //- Returns element in edgeIndices with minimum length
        static label minEdge
        (
            const triSurface& surf,
            const labelList& edgeIndices
        );

        //- Returns element in edgeIndices with minimum length
        static label maxEdge
        (
            const triSurface& surf,
            const labelList& edgeIndices
        );

        //- Merge points within distance
        static triSurface mergePoints
        (
            const triSurface& surf,
            const scalar mergeTol
        );

        //- Triangle (unit) normal. If nearest point to triangle on edge use
        //  edge normal (calculated on the fly); if on vertex use vertex normal.
        //  Uses planarTol.
        static vector surfaceNormal
        (
            const triSurface& surf,
446
            const label nearestFacei,
447 448 449
            const point& nearestPt
        );

450
        //- On which side of surface
451 452 453 454 455 456 457
        enum sideType
        {
            UNKNOWN,    // cannot be determined (e.g. non-manifold)
            INSIDE,     // inside
            OUTSIDE     // outside
        };

458
        //- If nearest point is on edgeI, determine on which side of surface
459 460 461 462 463 464 465 466 467 468 469
        //  sample is.
        static sideType edgeSide
        (
            const triSurface& surf,
            const point& sample,
            const point& nearestPoint,
            const label edgeI
        );

        //- Given nearest point (to sample) on surface determines which side
        //  sample is. Uses either face normal, edge normal or point normal
470
        //  (non-trivial). Uses triangle::classify.
471 472 473 474
        static sideType surfaceSide
        (
            const triSurface& surf,
            const point& sample,
475
            const label nearestFacei
476 477 478 479 480 481 482
        );

    // Triangulation of faces

        //- Simple triangulation of (selected patches of) boundaryMesh. Needs
        //  polyMesh (or polyBoundaryMesh) since only at this level are the
        //  triangles on neighbouring patches connected.
sergio's avatar
ENH:  
sergio committed
483
        //  Return faceMap from triI to faceI
484 485 486 487
        static triSurface triangulate
        (
            const polyBoundaryMesh& mBesh,
            const labelHashSet& includePatches,
sergio's avatar
ENH:  
sergio committed
488
            labelList& faceMap,
489 490 491
            const bool verbose = false
        );

492 493 494 495 496 497 498 499 500 501

        static triSurface triangulate
        (
            const polyBoundaryMesh& bMesh,
            const labelHashSet& includePatches,
            const boundBox& bBox,
            const bool verbose = false
        );


502 503 504 505 506 507 508 509 510 511 512 513 514 515
        //- Face-centre triangulation of (selected patches of) boundaryMesh.
        //  Needs
        //  polyMesh (or polyBoundaryMesh) since only at this level are the
        //  triangles on neighbouring patches connected.
        triSurface triangulateFaceCentre
        (
            const polyBoundaryMesh& mBesh,
            const labelHashSet& includePatches,
            const bool verbose = false
        );


    // Triangulation and interpolation

516 517 518 519
        //- Do unconstrained Delaunay of points. Returns triSurface with 3D
        //  points with z=0. All triangles in region 0.
        static triSurface delaunay2D(const List<vector2D>&);

520 521 522 523
        //- Calculate linear interpolation weights for point (guaranteed to be
        //  inside triangle)
        static void calcInterpolationWeights
        (
524 525
            const triPointRef& tri,
            const point& p,
526 527 528 529 530 531 532 533 534 535 536 537 538 539
            FixedList<scalar, 3>& weights
        );

        // Calculate weighting factors from samplePts to triangle it is in.
        // Uses linear search to find triangle.
        // Vertices are:
        //   (a b c)  : vertices of the triangle abc the point is in
        // or if the point is outside all triangles:
        //   (a b -1) : the edge ab the point is nearest to.
        //   (a -1 -1) : the vertex a the point is nearest to
        static void calcInterpolationWeights
        (
            const triSurface& s,
            const pointField& samplePts,
540 541
            List<FixedList<label, 3>>& allVerts,
            List<FixedList<scalar, 3>>& allWeights
542 543
        );

544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619

    // Curvature

        //- Weighting for normals of faces attached to vertex
        static scalar vertexNormalWeight
        (
            const triFace& f,
            const label pI,
            const vector& fN,
            const UList<point>& points
        );

        //- Weighted average of normals of attached faces
        static tmp<vectorField> vertexNormals(const triSurface& surf);

        //- Local coordinate-system for each point normal
        static tmp<triadField> vertexTriads
        (
            const triSurface& surf,
            const vectorField& pointNormals
        );

        //- Surface curvatures at the vertex points
        static tmp<scalarField> curvatures
        (
            const triSurface& surf,
            const vectorField& pointNormals,
            const triadField& pointTriads
        );

        //- Surface curvatures at the vertex points
        static tmp<scalarField> curvatures
        (
            const triSurface& surf
        );

        //- Write surface curvature at the vertex points and return the field
        static tmp<scalarField> writeCurvature
        (
            const Time& runTime,
            const word& basename,
            const triSurface& surf
        );


    // Closeness

        //- Check and write internal/external closeness fields
        static Pair<tmp<scalarField>> writeCloseness
        (
            const Time& runTime,
            const word& basename,
            const triSurface& surf,
            const scalar internalAngleTolerance = 45,
            const scalar externalAngleTolerance = 10
        );


    // Feature Proximity

        //- Calculate feature proximity
        scalarField featureProximity
        (
            const triSurface& surf,
            const scalar searchDistance
        );

        //- Check and write internal/external closeness fields
        static void writeFeatureProximity
        (
            const Time& runTime,
            const word& basename,
            const triSurface& surf,
            const bool writeVTK,
            const scalar searchDistance
        );
620 621


622 623 624
    // Surface checking functionality

        //- Check single triangle for (topological) validity
625 626 627 628 629 630
        static bool validTri
        (
            const triSurface&,
            const label facei,
            const bool verbose = true
        );
631 632

        //- Check single triangle for (topological) validity
633 634 635 636 637 638
        static bool validTri
        (
            const MeshedSurface<face>&,
            const label facei,
            const bool verbose = true
        );
639 640


641 642
    // Tracking

643
        //- Test point on plane of triangle to see if on edge or point or inside
644 645 646 647 648 649 650
        static surfaceLocation classify
        (
            const triSurface&,
            const label triI,
            const point& trianglePoint
        );

651 652
        //- Track on surface to get closer to point.
        //  Possible situations:
653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699
        //  - 1. reached endpoint
        //  - 2. reached edge (normal situation)
        //  - 3. reached end of surface (edge on single face)
        //  Input:
        //  - starting position+triangle/edge/point (so has to be on surface!)
        //  - (optional) previous position+triangle/edge/point to prevent
        //    going back. Set index (of triangle/edge/point) to -1 if not
        //    used.
        //  - end position+triangle/edge/point (so has to be on surface!)
        //  - plane to follow. Has to go through end point!
        //  Returns:
        //  - true if end point reached (situation 1)
        //  - new position+triangle/edge/point
        //  Caller has to check for situation 3 by checking that triangle()
        //  is not set.
        static surfaceLocation trackToEdge
        (
            const triSurface&,
            const surfaceLocation& start,
            const surfaceLocation& end,
            const plane& cutPlane
        );

        //- Track from edge to edge across surface. Uses trackToEdge.
        //  Not really useful by itself, more example of how to use trackToEdge.
        //  endInfo should be location on surface.
        //  hitInfo should be initialised to starting location (on surface as
        //  well). Upon return is set to end location.
        static void track
        (
            const triSurface&,
            const surfaceLocation& endInfo,
            const plane& cutPlane,
            surfaceLocation& hitInfo
        );
};


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

} // End namespace Foam

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#endif

// ************************************************************************* //