File: foamyHexMeshSurfaceSimplify.C

package info (click to toggle)
openfoam 4.1%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 163,028 kB
  • ctags: 58,990
  • sloc: cpp: 830,760; sh: 10,227; ansic: 8,215; xml: 745; lex: 437; awk: 194; sed: 91; makefile: 77; python: 18
file content (640 lines) | stat: -rw-r--r-- 17,730 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
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
137
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
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
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
261
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
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
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
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
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
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
/*---------------------------------------------------------------------------*\
 =========                   |
 \\      /   F ield          | OpenFOAM: The Open Source CFD Toolbox
  \\    /    O peration      |
   \\  /     A nd            | Copyright (C) 2012-2016 OpenFOAM Foundation
    \\/      M anipulation   |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    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.

    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
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Application
    foamyHexMeshSurfaceSimplify

Description
    Simplifies surfaces by resampling.

    Uses Thomas Lewiner's topology preserving MarchingCubes.
    (http://zeus.mat.puc-rio.br/tomlew/tomlew_uk.php)

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

#include "argList.H"
#include "Time.H"
#include "searchableSurfaces.H"
#include "conformationSurfaces.H"
#include "triSurfaceMesh.H"

#include "opt_octree.h"
#include "cube.h"

using namespace Foam;

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

class pointConversion
{
    const vector scale_;

    const vector offset_;

public:

    //- Construct from components
    pointConversion
    (
        const vector scale,
        const vector offset
    )
    :
        scale_(scale),
        offset_(offset)
    {}

    inline Point toLocal(const Foam::point& pt) const
    {
        Foam::point p = cmptMultiply(scale_, (pt + offset_));
        return Point(p.x(), p.y(), p.z());
    }

    inline Foam::point toGlobal(const Point& pt) const
    {
        point p(pt.x(), pt.y(), pt.z());
        return cmptDivide(p, scale_) - offset_;
    }
};




// For use in Fast-Dual Octree from Thomas Lewiner
class distanceCalc
:
    public ::data_access
{

    const Level min_level_;

    const conformationSurfaces& geometryToConformTo_;

    const pointConversion& converter_;


    // Private Member Functions

    scalar signedDistance(const Foam::point& pt) const
    {
        const searchableSurfaces& geometry = geometryToConformTo_.geometry();
        const labelList& surfaces = geometryToConformTo_.surfaces();

        static labelList nearestSurfaces;
        static scalarField distance;

        static pointField samples(1);
        samples[0] = pt;

        searchableSurfacesQueries::signedDistance
        (
            geometry,
            surfaces,
            samples,
            scalarField(1, GREAT),
            volumeType::OUTSIDE,
            nearestSurfaces,
            distance
        );

        return distance[0];
    }


public:

    // Constructors

        //- Construct from components
        distanceCalc
        (
            Level max_level_,
            real iso_val_,
            Level min_level,
            const conformationSurfaces& geometryToConformTo,
            const pointConversion& converter
        )
        :
            data_access(max_level_,iso_val_),
            min_level_(min_level),
            geometryToConformTo_(geometryToConformTo),
            converter_(converter)
        {}


    //- Destructor
    virtual ~distanceCalc()
    {}


    // Member Functions

        //- Test function
        virtual bool need_refine( const Cube &c )
        {
            int l = c.lv() ;

            if ( l >= _max_level ) return false;
            if ( l < min_level_ ) return true;

            treeBoundBox bb
            (
                converter_.toGlobal
                (
                    Point
                    (
                        c.xmin(),
                        c.ymin(),
                        c.zmin()
                    )
                ),
                converter_.toGlobal
                (
                    Point
                    (
                        c.xmax(),
                        c.ymax(),
                        c.zmax()
                    )
                )
            );

            const searchableSurfaces& geometry =
                geometryToConformTo_.geometry();
            const labelList& surfaces =
                geometryToConformTo_.surfaces();


            //- Uniform refinement around surface
            {
                forAll(surfaces, i)
                {
                    if (geometry[surfaces[i]].overlaps(bb))
                    {
                        return true;
                    }
                }
                return false;
            }


            ////- Surface intersects bb (but not using intersection test)
            //scalar ccDist = signedDistance(bb.midpoint());
            //scalar ccVal = ccDist - _iso_val;
            //if (mag(ccVal) < SMALL)
            //{
            //    return true;
            //}
            //const pointField points(bb.points());
            //forAll(points, pointi)
            //{
            //    scalar pointVal = signedDistance(points[pointi]) - _iso_val;
            //    if (ccVal*pointVal < 0)
            //    {
            //        return true;
            //    }
            //}
            //return false;


            ////- Refinement based on intersection with multiple planes.
            ////  Does not work well - too high a ratio between
            ////  neighbouring cubes.
            //const pointField points(bb.points());
            //const edgeList& edges = treeBoundBox::edges;
            //pointField start(edges.size());
            //pointField end(edges.size());
            //forAll(edges, i)
            //{
            //    start[i] = points[edges[i][0]];
            //    end[i] = points[edges[i][1]];
            //}
            //Foam::List<Foam::List<pointIndexHit>> hitInfo;
            //labelListList hitSurfaces;
            //searchableSurfacesQueries::findAllIntersections
            //(
            //    geometry,
            //    surfaces,
            //    start,
            //    end,
            //    hitSurfaces,
            //    hitInfo
            //);
            //
            //// Count number of intersections
            //label nInt = 0;
            //forAll(hitSurfaces, edgeI)
            //{
            //    nInt += hitSurfaces[edgeI].size();
            //}
            //
            //if (nInt == 0)
            //{
            //    // No surface intersected. See if there is one inside
            //    forAll(surfaces, i)
            //    {
            //        if (geometry[surfaces[i]].overlaps(bb))
            //        {
            //            return true;
            //        }
            //    }
            //    return false;
            //}
            //
            //// Check multiple surfaces
            //label baseSurfI = -1;
            //forAll(hitSurfaces, edgeI)
            //{
            //    const labelList& hSurfs = hitSurfaces[edgeI];
            //    forAll(hSurfs, i)
            //    {
            //        if (baseSurfI == -1)
            //        {
            //            baseSurfI = hSurfs[i];
            //        }
            //        else if (baseSurfI != hSurfs[i])
            //        {
            //            // Multiple surfaces
            //            return true;
            //        }
            //    }
            //}
            //
            //// Get normals
            //DynamicList<pointIndexHit> baseInfo(nInt);
            //forAll(hitInfo, edgeI)
            //{
            //    const Foam::List<pointIndexHit>& hits = hitInfo[edgeI];
            //    forAll(hits, i)
            //    {
            //        (void)hits[i].hitPoint();
            //        baseInfo.append(hits[i]);
            //    }
            //}
            //vectorField normals;
            //geometry[surfaces[baseSurfI]].getNormal(baseInfo, normals);
            //for (label i = 1; i < normals.size(); ++i)
            //{
            //    if ((normals[0] & normals[i]) < 0.9)
            //    {
            //        return true;
            //    }
            //}
            //labelList regions;
            //geometry[surfaces[baseSurfI]].getRegion(baseInfo, regions);
            //for (label i = 1; i < regions.size(); ++i)
            //{
            //    if (regions[0] != regions[i])
            //    {
            //        return true;
            //    }
            //}
            //return false;



            //samples[0] = point(c.xmin(), c.ymin(), c.zmin());
            //samples[1] = point(c.xmax(), c.ymin(), c.zmin());
            //samples[2] = point(c.xmax(), c.ymax(), c.zmin());
            //samples[3] = point(c.xmin(), c.ymax(), c.zmin());
            //
            //samples[4] = point(c.xmin(), c.ymin(), c.zmax());
            //samples[5] = point(c.xmax(), c.ymin(), c.zmax());
            //samples[6] = point(c.xmax(), c.ymax(), c.zmax());
            //samples[7] = point(c.xmin(), c.ymax(), c.zmax());

            //scalarField nearestDistSqr(8, GREAT);
            //
            //Foam::List<pointIndexHit> nearestInfo;
            //surf_.findNearest(samples, nearestDistSqr, nearestInfo);
            //vectorField normals;
            //surf_.getNormal(nearestInfo, normals);
            //
            //for (label i = 1; i < normals.size(); ++i)
            //{
            //    if ((normals[0] & normals[i]) < 0.5)
            //    {
            //        return true;
            //    }
            //}
            //return false;

            //// Check if surface octree same level
            //const labelList elems(surf_.tree().findBox(bb));
            //
            //if (elems.size() > 1)
            //{
            //    return true;
            //}
            //else
            //{
            //  return false;
            //}
        }

        //- Data function
        virtual real value_at( const Cube &c )
        {
            return signedDistance(converter_.toGlobal(c)) - _iso_val;
        }
};


// Main program:

int main(int argc, char *argv[])
{
    argList::addNote
    (
        "Re-sample surfaces used in foamyHexMesh operation"
    );
    argList::validArgs.append("outputName");

    #include "setRootCase.H"
    #include "createTime.H"
    runTime.functionObjects().off();

    const fileName exportName = args.args()[1];

    Info<< "Reading surfaces as specified in the foamyHexMeshDict and"
        << " writing a re-sampled surface to " << exportName
        << nl << endl;

    cpuTime timer;

    IOdictionary foamyHexMeshDict
    (
        IOobject
        (
            "foamyHexMeshDict",
            runTime.system(),
            runTime,
            IOobject::MUST_READ_IF_MODIFIED,
            IOobject::NO_WRITE
        )
    );

    // Define/load all geometry
    searchableSurfaces allGeometry
    (
        IOobject
        (
            "cvSearchableSurfaces",
            runTime.constant(),
            "triSurface",
            runTime,
            IOobject::MUST_READ,
            IOobject::NO_WRITE
        ),
        foamyHexMeshDict.subDict("geometry"),
        foamyHexMeshDict.lookupOrDefault("singleRegionName", true)
    );

    Info<< "Geometry read in = "
        << timer.cpuTimeIncrement() << " s." << nl << endl;


    Random rndGen(64293*Pstream::myProcNo());

    conformationSurfaces geometryToConformTo
    (
        runTime,
        rndGen,
        allGeometry,
        foamyHexMeshDict.subDict("surfaceConformation")
    );

    Info<< "Set up geometry in = "
        << timer.cpuTimeIncrement() << " s." << nl << endl;


    const searchableSurfaces& geometry = geometryToConformTo.geometry();
    const labelList& surfaces = geometryToConformTo.surfaces();


    const label minLevel = 2;

    // The max cube size follows from the minLevel and the default cube size
    // (1)
    const scalar maxSize = 1.0 / (1 << minLevel);
    const scalar halfMaxSize = 0.5*maxSize;


    // Scale the geometry to fit within
    // halfMaxSize .. 1-halfMaxSize

    scalar wantedRange = 1.0-maxSize;

    const treeBoundBox bb = geometryToConformTo.globalBounds();

    const vector scale = cmptDivide
    (
        point(wantedRange, wantedRange, wantedRange),
        bb.span()
    );
    const vector offset =
        cmptDivide
        (
            point(halfMaxSize, halfMaxSize, halfMaxSize),
            scale
        )
       -bb.min();


    const pointConversion converter(scale, offset);


    // Marching cubes

    OptOctree octree;

    distanceCalc ref
    (
        8,          //maxLevel
        0.0,        //distance
        minLevel,   //minLevel
        geometryToConformTo,
        converter
    );

    octree.refine(&ref);
    octree.set_impl(&ref);

    Info<< "Calculated octree in = "
        << timer.cpuTimeIncrement() << " s." << nl << endl;

    MarchingCubes& mc = octree.mc();

    mc.clean_all() ;
    octree.build_isosurface(&ref) ;

    Info<< "Constructed iso surface of distance in = "
        << timer.cpuTimeIncrement() << " s." << nl << endl;

    // Write output file
    if (mc.ntrigs() > 0)
    {
        Triangle* triangles = mc.triangles();
        label nTris = mc.ntrigs();
        Foam::DynamicList<labelledTri> tris(mc.ntrigs());
        for (label triI = 0; triI < nTris; ++triI)
        {
            const Triangle& t = triangles[triI];
            if (t.v1 != t.v2 && t.v1 != t.v3 && t.v2 != t.v3)
            {
                tris.append
                (
                    labelledTri
                    (
                        triangles[triI].v1,
                        triangles[triI].v2,
                        triangles[triI].v3,
                        0                       // region
                    )
                );
            }
        }


        Point* vertices = mc.vertices();
        pointField points(mc.nverts());
        forAll(points, pointi)
        {
            const Point& v = vertices[pointi];
            points[pointi] = converter.toGlobal(v);
        }


        // Find correspondence to original surfaces
        labelList regionOffsets(surfaces.size());
        label nRegions = 0;
        forAll(surfaces, i)
        {
            const wordList& regions = geometry[surfaces[i]].regions();
            regionOffsets[i] = nRegions;
            nRegions += regions.size();
        }


        geometricSurfacePatchList patches(nRegions);
        nRegions = 0;
        forAll(surfaces, i)
        {
            const wordList& regions = geometry[surfaces[i]].regions();

            forAll(regions, regionI)
            {
                patches[nRegions] = geometricSurfacePatch
                (
                    "patch",
                    geometry[surfaces[i]].name() + "_" + regions[regionI],
                    nRegions
                );
                nRegions++;
            }
        }

        triSurface s(tris, patches, points, true);
        tris.clearStorage();

        Info<< "Extracted triSurface in = "
            << timer.cpuTimeIncrement() << " s." << nl << endl;


        // Find out region on local surface of nearest point
        {
            Foam::List<pointIndexHit> hitInfo;
            labelList hitSurfaces;
            geometryToConformTo.findSurfaceNearest
            (
                s.faceCentres(),
                scalarField(s.size(), sqr(GREAT)),
                hitInfo,
                hitSurfaces
            );

            // Get region
            DynamicList<pointIndexHit> surfInfo(hitSurfaces.size());
            DynamicList<label> surfIndices(hitSurfaces.size());

            forAll(surfaces, surfI)
            {
                // Extract info on this surface
                surfInfo.clear();
                surfIndices.clear();
                forAll(hitSurfaces, triI)
                {
                    if (hitSurfaces[triI] == surfI)
                    {
                        surfInfo.append(hitInfo[triI]);
                        surfIndices.append(triI);
                    }
                }

                // Calculate sideness of these surface points
                labelList region;
                geometry[surfaces[surfI]].getRegion(surfInfo, region);

                forAll(region, i)
                {
                    label triI = surfIndices[i];
                    s[triI].region() = regionOffsets[surfI]+region[i];
                }
            }
        }

        Info<< "Re-patched surface in = "
            << timer.cpuTimeIncrement() << " s." << nl << endl;

        triSurfaceMesh smesh
        (
            IOobject
            (
                exportName,
                runTime.constant(), // instance
                "triSurface",
                runTime,            // registry
                IOobject::NO_READ,
                IOobject::AUTO_WRITE,
                false
            ),
            s
        );

        Info<< "writing surfMesh:\n  "
            << smesh.searchableSurface::objectPath() << nl << endl;
        smesh.searchableSurface::write();

        Info<< "Written surface in = "
            << timer.cpuTimeIncrement() << " s." << nl << endl;
    }

    mc.clean_all() ;

    Info<< "End\n" << endl;

    return 0;
}


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