surfaceLambdaMuSmooth.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 {
131  #include "removeCaseOptions.H"
132 
133  argList::validOptions.clear();
134  argList::validArgs.append("surface file");
135  argList::validArgs.append("output surface file");
136  argList::validArgs.append("lambda (0..1)");
137  argList::validArgs.append("mu (0..1)");
138  argList::validArgs.append("iterations");
140  (
141  "featureFile",
142  "fix points from a file containing feature points and edges"
143  );
144  argList args(argc, argv);
145 
146  const fileName surfFileName = args[1];
147  const fileName outFileName = args[2];
148  const scalar lambda = args.argRead<scalar>(3);
149  const scalar mu = args.argRead<scalar>(4);
150  const label iters = args.argRead<label>(5);
151 
152  if (lambda < 0 || lambda > 1)
153  {
155  << lambda << endl
156  << "0: no change 1: move vertices to average of neighbours"
157  << exit(FatalError);
158  }
159  if (mu < 0 || mu > 1)
160  {
162  << mu << endl
163  << "0: no change 1: move vertices to average of neighbours"
164  << exit(FatalError);
165  }
166 
167  Info<< "lambda : " << lambda << nl
168  << "mu : " << mu << nl
169  << "Iters : " << iters << nl
170  << "Reading surface from " << surfFileName << " ..." << endl;
171 
172  meshedSurface surf1(surfFileName);
173 
174  Info<< "Faces : " << surf1.size() << nl
175  << "Vertices : " << surf1.nPoints() << nl
176  << "Bounding Box: " << boundBox(surf1.localPoints()) << endl;
177 
178  PackedBoolList fixedPoints(surf1.localPoints().size(), false);
179 
180  if (args.optionFound("featureFile"))
181  {
182  const fileName featureFileName(args.option("featureFile"));
183  Info<< "Reading features from " << featureFileName << " ..." << endl;
184 
185  edgeMesh feMesh(featureFileName);
186 
187  getFixedPoints(feMesh, surf1.localPoints(), fixedPoints);
188 
189  Info<< "Number of fixed points on surface = " << fixedPoints.count()
190  << endl;
191  }
192 
193  pointField newPoints(surf1.localPoints());
194 
195  for (label iter = 0; iter < iters; iter++)
196  {
197  // Lambda
198  {
199  pointField newLocalPoints
200  (
201  (1 - lambda)*surf1.localPoints()
202  + lambda*avg(surf1, fixedPoints)
203  );
204 
205  pointField newPoints(surf1.points());
206  UIndirectList<point>(newPoints, surf1.meshPoints()) =
207  newLocalPoints;
208 
209  surf1.movePoints(newPoints);
210  }
211 
212  // Mu
213  if (mu != 0)
214  {
215  pointField newLocalPoints
216  (
217  (1 + mu)*surf1.localPoints()
218  - mu*avg(surf1, fixedPoints)
219  );
220 
221  pointField newPoints(surf1.points());
222  UIndirectList<point>(newPoints, surf1.meshPoints()) =
223  newLocalPoints;
224 
225  surf1.movePoints(newPoints);
226  }
227  }
228 
229  Info<< "Writing surface to " << outFileName << " ..." << endl;
230  surf1.write(outFileName);
231 
232  Info<< "End\n" << endl;
233 
234  return 0;
235 }
236 
237 
238 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual Ostream & write(const char)=0
Write character.
A bit-packed bool list.
unsigned int count() const
Count number of bits set, O(log(n))
Definition: PackedList.C:55
A List with indirect addressing.
Definition: UIndirectList.H:60
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:103
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
const string & option(const word &opt) const
Return the argument string associated with the named option.
Definition: argListI.H:108
static HashTable< string > validOptions
A list of valid options.
Definition: argList.H:156
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:183
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
Points connected by edges.
Definition: edgeMesh.H:72
const pointField & points() const
Return points.
Definition: edgeMeshI.H:62
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A class for handling file names.
Definition: fileName.H:82
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar lambda(viscosity->lookup("lambda"))
Determine correspondence between points. See below.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar mu
Atomic mass unit.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:105
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
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
error FatalError
static const char nl
Definition: Ostream.H:260
Foam::argList args(argc, argv)