surfaceLambdaMuSmooth.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  surfaceLambdaMuSmooth
26 
27 Description
28  Smooths a surface using lambda/mu smoothing.
29 
30  To get laplacian smoothing, set lambda to the relaxation factor and mu to
31  zero.
32 
33  Provide an edgeMesh file containing points that are not to be moved during
34  smoothing in order to preserve features.
35 
36  lambda/mu smoothing: G. Taubin, IBM Research report Rc-19923 (02/01/95)
37  "A signal processing approach to fair surface design"
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "argList.H"
42 #include "boundBox.H"
43 #include "edgeMesh.H"
44 #include "matchPoints.H"
45 #include "MeshedSurfaces.H"
46 
47 using namespace Foam;
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
52 (
53  const meshedSurface& s,
54  const PackedBoolList& fixedPoints
55 )
56 {
57  const labelListList& pointEdges = s.pointEdges();
58 
59  tmp<pointField> tavg(new pointField(s.nPoints(), Zero));
60  pointField& avg = tavg.ref();
61 
62  forAll(pointEdges, vertI)
63  {
64  vector& avgPos = avg[vertI];
65 
66  if (fixedPoints[vertI])
67  {
68  avgPos = s.localPoints()[vertI];
69  }
70  else
71  {
72  const labelList& pEdges = pointEdges[vertI];
73 
74  forAll(pEdges, myEdgeI)
75  {
76  const edge& e = s.edges()[pEdges[myEdgeI]];
77 
78  label otherVertI = e.otherVertex(vertI);
79 
80  avgPos += s.localPoints()[otherVertI];
81  }
82 
83  avgPos /= pEdges.size();
84  }
85  }
86 
87  return tavg;
88 }
89 
90 
91 void getFixedPoints
92 (
93  const edgeMesh& feMesh,
94  const pointField& points,
95  PackedBoolList& fixedPoints
96 )
97 {
98  scalarList matchDistance(feMesh.points().size(), 1e-1);
99  labelList from0To1;
100 
101  bool matchedAll = matchPoints
102  (
103  feMesh.points(),
104  points,
105  matchDistance,
106  false,
107  from0To1
108  );
109 
110  if (!matchedAll)
111  {
113  << "Did not match all feature points to points on the surface"
114  << endl;
115  }
116 
117  forAll(from0To1, fpI)
118  {
119  if (from0To1[fpI] != -1)
120  {
121  fixedPoints[from0To1[fpI]] = true;
122  }
123  }
124 }
125 
126 
127 // Main program:
128 
129 int main(int argc, char *argv[])
130 {
132  argList::validOptions.clear();
133  argList::validArgs.append("surfaceFile");
134  argList::validArgs.append("lambda (0..1)");
135  argList::validArgs.append("mu (0..1)");
136  argList::validArgs.append("iterations");
137  argList::validArgs.append("output surfaceFile");
139  (
140  "featureFile",
141  "fix points from a file containing feature points and edges"
142  );
143  argList args(argc, argv);
144 
145  const fileName surfFileName = args[1];
146  const scalar lambda = args.argRead<scalar>(2);
147  const scalar mu = args.argRead<scalar>(3);
148  const label iters = args.argRead<label>(4);
149  const fileName outFileName = args[5];
150 
151  if (lambda < 0 || lambda > 1)
152  {
154  << lambda << endl
155  << "0: no change 1: move vertices to average of neighbours"
156  << exit(FatalError);
157  }
158  if (mu < 0 || mu > 1)
159  {
161  << mu << endl
162  << "0: no change 1: move vertices to average of neighbours"
163  << exit(FatalError);
164  }
165 
166  Info<< "lambda : " << lambda << nl
167  << "mu : " << mu << nl
168  << "Iters : " << iters << nl
169  << "Reading surface from " << surfFileName << " ..." << endl;
170 
171  meshedSurface surf1(surfFileName);
172 
173  Info<< "Faces : " << surf1.size() << nl
174  << "Vertices : " << surf1.nPoints() << nl
175  << "Bounding Box: " << boundBox(surf1.localPoints()) << endl;
176 
177  PackedBoolList fixedPoints(surf1.localPoints().size(), false);
178 
179  if (args.optionFound("featureFile"))
180  {
181  const fileName featureFileName(args.option("featureFile"));
182  Info<< "Reading features from " << featureFileName << " ..." << endl;
183 
184  edgeMesh feMesh(featureFileName);
185 
186  getFixedPoints(feMesh, surf1.localPoints(), fixedPoints);
187 
188  Info<< "Number of fixed points on surface = " << fixedPoints.count()
189  << endl;
190  }
191 
192  pointField newPoints(surf1.localPoints());
193 
194  for (label iter = 0; iter < iters; iter++)
195  {
196  // Lambda
197  {
198  pointField newLocalPoints
199  (
200  (1 - lambda)*surf1.localPoints()
201  + lambda*avg(surf1, fixedPoints)
202  );
203 
204  pointField newPoints(surf1.points());
205  UIndirectList<point>(newPoints, surf1.meshPoints()) =
206  newLocalPoints;
207 
208  surf1.movePoints(newPoints);
209  }
210 
211  // Mu
212  if (mu != 0)
213  {
214  pointField newLocalPoints
215  (
216  (1 + mu)*surf1.localPoints()
217  - mu*avg(surf1, fixedPoints)
218  );
219 
220  pointField newPoints(surf1.points());
221  UIndirectList<point>(newPoints, surf1.meshPoints()) =
222  newLocalPoints;
223 
224  surf1.movePoints(newPoints);
225  }
226  }
227 
228  Info<< "Writing surface to " << outFileName << " ..." << endl;
229  surf1.write(outFileName);
230 
231  Info<< "End\n" << endl;
232 
233  return 0;
234 }
235 
236 
237 // ************************************************************************* //
unsigned int count() const
Count number of bits set, O(log(n))
Definition: PackedList.C:55
label nPoints() const
Return number of points supporting patch faces.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
label otherVertex(const label a) const
Given one vertex, return the other.
Definition: edgeI.H:103
Determine correspondence between points. See below.
const labelListList & pointEdges() const
Return point-edge addressing.
const pointField & points() const
Return points.
Definition: edgeMeshI.H:39
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition: matchPoints.C:33
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const pointField & points
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:93
static const zero Zero
Definition: zero.H:91
static const char nl
Definition: Ostream.H:262
Points connected by edges.
Definition: edgeMesh.H:69
const string & option(const word &opt) const
Return the argument string associated with the named option.
Definition: argListI.H:102
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:177
A bit-packed bool list.
#define WarningInFunction
Report a warning using Foam::Warning.
A List with indirect addressing.
Definition: fvMatrix.H:106
messageStream Info
virtual Ostream & write(const token &)=0
Write next token to stream.
A class for managing temporary objects.
Definition: PtrList.H:54
Foam::argList args(argc, argv)
static HashTable< string > validOptions
A list of valid options.
Definition: argList.H:157
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Namespace for OpenFOAM.