File: surfaceLambdaMuSmooth.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 (237 lines) | stat: -rw-r--r-- 6,455 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
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-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
    surfaceLambdaMuSmooth

Description
    Smooths a surface using lambda/mu smoothing.

    To get laplacian smoothing, set lambda to the relaxation factor and mu to
    zero.

    Provide an edgeMesh file containing points that are not to be moved during
    smoothing in order to preserve features.

    lambda/mu smoothing: G. Taubin, IBM Research report Rc-19923 (02/01/95)
    "A signal processing approach to fair surface design"

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

#include "argList.H"
#include "boundBox.H"
#include "edgeMesh.H"
#include "matchPoints.H"
#include "MeshedSurfaces.H"

using namespace Foam;

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

tmp<pointField> avg
(
    const meshedSurface& s,
    const PackedBoolList& fixedPoints
)
{
    const labelListList& pointEdges = s.pointEdges();

    tmp<pointField> tavg(new pointField(s.nPoints(), Zero));
    pointField& avg = tavg.ref();

    forAll(pointEdges, vertI)
    {
        vector& avgPos = avg[vertI];

        if (fixedPoints[vertI])
        {
            avgPos = s.localPoints()[vertI];
        }
        else
        {
            const labelList& pEdges = pointEdges[vertI];

            forAll(pEdges, myEdgeI)
            {
                const edge& e = s.edges()[pEdges[myEdgeI]];

                label otherVertI = e.otherVertex(vertI);

                avgPos += s.localPoints()[otherVertI];
            }

            avgPos /= pEdges.size();
        }
    }

    return tavg;
}


void getFixedPoints
(
    const edgeMesh& feMesh,
    const pointField& points,
    PackedBoolList& fixedPoints
)
{
    scalarList matchDistance(feMesh.points().size(), 1e-1);
    labelList from0To1;

    bool matchedAll = matchPoints
    (
        feMesh.points(),
        points,
        matchDistance,
        false,
        from0To1
    );

    if (!matchedAll)
    {
        WarningInFunction
            << "Did not match all feature points to points on the surface"
            << endl;
    }

    forAll(from0To1, fpI)
    {
        if (from0To1[fpI] != -1)
        {
            fixedPoints[from0To1[fpI]] = true;
        }
    }
}


// Main program:

int main(int argc, char *argv[])
{
    argList::noParallel();
    argList::validOptions.clear();
    argList::validArgs.append("surfaceFile");
    argList::validArgs.append("lambda (0..1)");
    argList::validArgs.append("mu (0..1)");
    argList::validArgs.append("iterations");
    argList::validArgs.append("output surfaceFile");
    argList::addOption
    (
        "featureFile",
        "fix points from a file containing feature points and edges"
    );
    argList args(argc, argv);

    const fileName surfFileName = args[1];
    const scalar lambda = args.argRead<scalar>(2);
    const scalar mu = args.argRead<scalar>(3);
    const label  iters = args.argRead<label>(4);
    const fileName outFileName = args[5];

    if (lambda < 0 || lambda > 1)
    {
        FatalErrorInFunction
            << lambda << endl
            << "0: no change   1: move vertices to average of neighbours"
            << exit(FatalError);
    }
    if (mu < 0 || mu > 1)
    {
        FatalErrorInFunction
            << mu << endl
            << "0: no change   1: move vertices to average of neighbours"
            << exit(FatalError);
    }

    Info<< "lambda      : " << lambda << nl
        << "mu          : " << mu << nl
        << "Iters       : " << iters << nl
        << "Reading surface from " << surfFileName << " ..." << endl;

    meshedSurface surf1(surfFileName);

    Info<< "Faces       : " << surf1.size() << nl
        << "Vertices    : " << surf1.nPoints() << nl
        << "Bounding Box: " << boundBox(surf1.localPoints()) << endl;

    PackedBoolList fixedPoints(surf1.localPoints().size(), false);

    if (args.optionFound("featureFile"))
    {
        const fileName featureFileName(args.option("featureFile"));
        Info<< "Reading features from " << featureFileName << " ..." << endl;

        edgeMesh feMesh(featureFileName);

        getFixedPoints(feMesh, surf1.localPoints(), fixedPoints);

        Info<< "Number of fixed points on surface = " << fixedPoints.count()
            << endl;
    }

    pointField newPoints(surf1.localPoints());

    for (label iter = 0; iter < iters; iter++)
    {
        // Lambda
        {
            pointField newLocalPoints
            (
                (1 - lambda)*surf1.localPoints()
              + lambda*avg(surf1, fixedPoints)
            );

            pointField newPoints(surf1.points());
            UIndirectList<point>(newPoints, surf1.meshPoints()) =
                newLocalPoints;

            surf1.movePoints(newPoints);
        }

        // Mu
        if (mu != 0)
        {
            pointField newLocalPoints
            (
                (1 + mu)*surf1.localPoints()
              - mu*avg(surf1, fixedPoints)
            );

            pointField newPoints(surf1.points());
            UIndirectList<point>(newPoints, surf1.meshPoints()) =
                newLocalPoints;

            surf1.movePoints(newPoints);
        }
    }

    Info<< "Writing surface to " << outFileName << " ..." << endl;
    surf1.write(outFileName);

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

    return 0;
}


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