triSurfaceCurvature.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "triSurfaceTools.H"
29
30#include "triSurface.H"
31#include "MeshedSurfaces.H"
32#include "triSurfaceFields.H"
33#include "OFstream.H"
34#include "plane.H"
35#include "tensor2D.H"
36#include "symmTensor2D.H"
37#include "scalarMatrices.H"
38#include "transform.H"
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
43(
44 const triFace& f,
45 const label pI,
46 const vector& fN,
47 const UList<point>& points
48)
49{
50 label index = f.find(pI);
51
52 if (index == -1)
53 {
55 << "Point not in face" << abort(FatalError);
56 }
57
58 const vector e1 = points[f[index]] - points[f[f.fcIndex(index)]];
59 const vector e2 = points[f[index]] - points[f[f.rcIndex(index)]];
60
61 return mag(fN)/(magSqr(e1)*magSqr(e2) + VSMALL);
62}
63
64
67{
68 // Weighted average of normals of faces attached to the vertex
69 // Weight = fA / (mag(e0)^2 * mag(e1)^2);
70
71 Info<< "Calculating vertex normals" << endl;
72
73 auto tpointNormals = tmp<vectorField>::New(surf.nPoints(), Zero);
74 auto& pointNormals = tpointNormals.ref();
75
76 const pointField& points = surf.points();
77 const labelListList& pointFaces = surf.pointFaces();
78 const labelList& meshPoints = surf.meshPoints();
79
80 forAll(pointFaces, pI)
81 {
82 const labelList& pFaces = pointFaces[pI];
83
84 for (const label facei : pFaces)
85 {
86 const triFace& f = surf[facei];
87
88 const vector areaNorm = f.areaNormal(points);
89
90 const scalar weight = vertexNormalWeight
91 (
92 f,
93 meshPoints[pI],
94 areaNorm,
95 points
96 );
97
98 pointNormals[pI] += weight * areaNorm;
99 }
100
101 pointNormals[pI].normalise();
102 }
103
104 return tpointNormals;
105}
106
107
110(
111 const triSurface& surf,
112 const vectorField& pointNormals
113)
114{
115 const pointField& points = surf.points();
116 const Map<label>& meshPointMap = surf.meshPointMap();
117
118 auto tpointTriads = tmp<triadField>::New(points.size());
119 auto& pointTriads = tpointTriads.ref();
120
121 forAll(points, pI)
122 {
123 const point& pt = points[pI];
124 const vector& normal = pointNormals[meshPointMap[pI]];
125
126 if (mag(normal) < SMALL)
127 {
128 pointTriads[meshPointMap[pI]] = triad::unset;
129 continue;
130 }
131
132 plane p(pt, normal);
133
134 // Pick arbitrary point in plane
135 vector dir1 = normalised(pt - p.somePointInPlane(1e-3));
136 vector dir2 = normalised(dir1 ^ normal);
137
138 pointTriads[meshPointMap[pI]] = triad(dir1, dir2, normal);
139 }
140
141 return tpointTriads;
142}
143
144
147(
148 const triSurface& surf,
149 const vectorField& pointNormals,
150 const triadField& pointTriads
151)
152{
153 Info<< "Calculating face curvature" << endl;
154
155 const pointField& points = surf.points();
156 const labelList& meshPoints = surf.meshPoints();
157 const Map<label>& meshPointMap = surf.meshPointMap();
158
159 List<symmTensor2D> pointFundamentalTensors
160 (
161 points.size(),
163 );
164
165 scalarList accumulatedWeights(points.size(), Zero);
166
167 forAll(surf, fI)
168 {
169 const triFace& f = surf[fI];
170 const edgeList fEdges = f.edges();
171
172 // Calculate the edge vectors and the normal differences
173 vectorField edgeVectors(f.size(), Zero);
174 vectorField normalDifferences(f.size(), Zero);
175
176 forAll(fEdges, feI)
177 {
178 const edge& e = fEdges[feI];
179
180 edgeVectors[feI] = e.vec(points);
181 normalDifferences[feI] =
182 pointNormals[meshPointMap[e[0]]]
183 - pointNormals[meshPointMap[e[1]]];
184 }
185
186 // Set up a local coordinate system for the face
187 const vector& e0 = edgeVectors[0];
188 const vector eN = f.areaNormal(points);
189 const vector e1 = (e0 ^ eN);
190
191 if (magSqr(eN) < ROOTVSMALL)
192 {
193 continue;
194 }
195
196 triad faceCoordSys(e0, e1, eN);
197 faceCoordSys.normalize();
198
199 // Construct the matrix to solve
201 scalarDiagonalMatrix Z(3, 0);
202
203 // Least Squares
204 for (label i = 0; i < 3; ++i)
205 {
206 scalar x = edgeVectors[i] & faceCoordSys[0];
207 scalar y = edgeVectors[i] & faceCoordSys[1];
208
209 T(0, 0) += sqr(x);
210 T(1, 0) += x*y;
211 T(1, 1) += sqr(x) + sqr(y);
212 T(2, 1) += x*y;
213 T(2, 2) += sqr(y);
214
215 scalar dndx = normalDifferences[i] & faceCoordSys[0];
216 scalar dndy = normalDifferences[i] & faceCoordSys[1];
217
218 Z[0] += dndx*x;
219 Z[1] += dndx*y + dndy*x;
220 Z[2] += dndy*y;
221 }
222
223 // Perform Cholesky decomposition and back substitution.
224 // Decomposed matrix is in T and solution is in Z.
225 LUsolve(T, Z);
226 symmTensor2D secondFundamentalTensor(Z[0], Z[1], Z[2]);
227
228 // Loop over the face points adding the contribution of the face
229 // curvature to the points.
230 forAll(f, fpI)
231 {
232 const label patchPointIndex = meshPointMap[f[fpI]];
233
234 const triad& ptCoordSys = pointTriads[patchPointIndex];
235
236 if (!ptCoordSys.set())
237 {
238 continue;
239 }
240
241 // Rotate faceCoordSys to ptCoordSys
242 tensor rotTensor = rotationTensor(ptCoordSys[2], faceCoordSys[2]);
243 triad rotatedFaceCoordSys = rotTensor & tensor(faceCoordSys);
244
245 // Project the face curvature onto the point plane
246
247 vector2D cmp1
248 (
249 ptCoordSys[0] & rotatedFaceCoordSys[0],
250 ptCoordSys[0] & rotatedFaceCoordSys[1]
251 );
252 vector2D cmp2
253 (
254 ptCoordSys[1] & rotatedFaceCoordSys[0],
255 ptCoordSys[1] & rotatedFaceCoordSys[1]
256 );
257
258 tensor2D projTensor
259 (
260 cmp1,
261 cmp2
262 );
263
264 symmTensor2D projectedFundamentalTensor
265 (
266 projTensor.x() & (secondFundamentalTensor & projTensor.x()),
267 projTensor.x() & (secondFundamentalTensor & projTensor.y()),
268 projTensor.y() & (secondFundamentalTensor & projTensor.y())
269 );
270
271 // Calculate weight
272 // TODO: Voronoi area weighting
274 (
275 f,
276 meshPoints[patchPointIndex],
277 f.areaNormal(points),
278 points
279 );
280
281 // Sum contribution of face to this point
282 pointFundamentalTensors[patchPointIndex] +=
283 weight*projectedFundamentalTensor;
284
285 accumulatedWeights[patchPointIndex] += weight;
286 }
287
288 if (false)
289 {
290 Info<< "Points = " << points[f[0]] << " "
291 << points[f[1]] << " "
292 << points[f[2]] << endl;
293 Info<< "edgeVecs = " << edgeVectors[0] << " "
294 << edgeVectors[1] << " "
295 << edgeVectors[2] << endl;
296 Info<< "normDiff = " << normalDifferences[0] << " "
297 << normalDifferences[1] << " "
298 << normalDifferences[2] << endl;
299 Info<< "faceCoordSys = " << faceCoordSys << endl;
300 Info<< "T = " << T << endl;
301 Info<< "Z = " << Z << endl;
302 }
303 }
304
305 auto tcurvatureAtPoints = tmp<scalarField>::New(points.size(), Zero);
306 scalarField& curvatureAtPoints = tcurvatureAtPoints.ref();
307
308 forAll(curvatureAtPoints, pI)
309 {
310 pointFundamentalTensors[pI] /= (accumulatedWeights[pI] + SMALL);
311
312 vector2D principalCurvatures(eigenValues(pointFundamentalTensors[pI]));
313
314 //scalar curvature =
315 // (principalCurvatures[0] + principalCurvatures[1])/2;
316 scalar curvature = max
317 (
318 mag(principalCurvatures[0]),
319 mag(principalCurvatures[1])
320 );
321 //scalar curvature = principalCurvatures[0]*principalCurvatures[1];
322
323 curvatureAtPoints[meshPoints[pI]] = curvature;
324 }
325
326 return tcurvatureAtPoints;
327}
328
329
332(
333 const triSurface& surf
334)
335{
337 tmp<triadField> triads = triSurfaceTools::vertexTriads(surf, norms());
338
339 tmp<scalarField> curv = curvatures(surf, norms(), triads());
340 norms.clear();
341 triads.clear();
342
343 return curv;
344}
345
346
349(
350 const Time& runTime,
351 const word& basename,
352 const triSurface& surf
353)
354{
355 Info<< "Extracting curvature of surface at the points." << endl;
356
358 scalarField& curv = tcurv.ref();
359
361 (
363 (
364 basename + ".curvature",
366 "triSurface",
367 runTime,
370 ),
371 surf,
372 dimLength,
374 );
375
376 outputField.swap(curv);
377 outputField.write();
378 outputField.swap(curv);
379
380 return tcurv;
381}
382
383
384// ************************************************************************* //
scalar y
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
label nPoints() const
Number of points supporting patch faces.
const Map< label > & meshPointMap() const
Mesh point map.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const labelListList & pointFaces() const
Return point-face addressing.
A templated (2 x 2) symmetric tensor of objects of <T>, effectively containing 3 elements,...
Definition: SymmTensor2D.H:61
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
A templated (2 x 2) tensor of objects of <T> derived from VectorSpace.
Definition: Tensor2D.H:59
Vector2D< Cmpt > y() const
Extract vector for row 1.
Definition: Tensor2DI.H:154
Vector2D< Cmpt > x() const
Extract vector for row 0.
Definition: Tensor2DI.H:148
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void swap(UList< T > &list)
Swap content with another UList of the same type in constant time.
Definition: UListI.H:434
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
label rcIndex(const label i) const noexcept
Definition: UListI.H:67
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition: VectorI.H:123
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:95
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:72
static tmp< scalarField > curvatures(const triSurface &surf, const vectorField &pointNormals, const triadField &pointTriads)
Surface curvatures at the vertex points.
static tmp< triadField > vertexTriads(const triSurface &surf, const vectorField &pointNormals)
Local coordinate-system for each point normal.
static scalar vertexNormalWeight(const triFace &f, const label pI, const vector &fN, const UList< point > &points)
Weighting for normals of faces attached to vertex.
static tmp< scalarField > writeCurvature(const Time &runTime, const word &basename, const triSurface &surf)
Write surface curvature at the vertex points and return the field.
static tmp< vectorField > vertexNormals(const triSurface &surf)
Weighted average of normals of attached faces.
Triangulated surface description with patch information.
Definition: triSurface.H:79
Representation of a 3D Cartesian coordinate system as a Vector of row vectors.
Definition: triad.H:67
static const triad unset
Definition: triad.H:97
bool set(const direction d) const
Is the vector in the direction d set.
Definition: triadI.H:70
void normalize()
Same as normalise.
Definition: triad.H:170
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedVector eigenValues(const dimensionedSymmTensor &dt)
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Tensor< scalar > tensor
Definition: symmTensor.H:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59
3D tensor transformation operations.
Fields for triSurface.