Loading...
Searching...
No Matches
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 const label index = f.find(pI);
51
52 if (index == -1)
53 {
55 << "Point not in face" << abort(FatalError);
56
57 return 0; // make gcc-15 happy by avoiding triFace[-1]
58 }
59
60 const vector e1 = points[f[index]] - points[f[f.fcIndex(index)]];
61 const vector e2 = points[f[index]] - points[f[f.rcIndex(index)]];
62
63 return mag(fN)/(magSqr(e1)*magSqr(e2) + VSMALL);
64}
65
66
69{
70 // Weighted average of normals of faces attached to the vertex
71 // Weight = fA / (mag(e0)^2 * mag(e1)^2);
72
73 Info<< "Calculating vertex normals" << endl;
74
75 auto tpointNormals = tmp<vectorField>::New(surf.nPoints(), Zero);
76 auto& pointNormals = tpointNormals.ref();
77
78 const pointField& points = surf.points();
79 const labelListList& pointFaces = surf.pointFaces();
80 const labelList& meshPoints = surf.meshPoints();
81
82 forAll(pointFaces, pI)
83 {
84 const labelList& pFaces = pointFaces[pI];
85
86 for (const label facei : pFaces)
87 {
88 const triFace& f = surf[facei];
89
90 const vector areaNorm = f.areaNormal(points);
91
92 const scalar weight = vertexNormalWeight
93 (
94 f,
95 meshPoints[pI],
96 areaNorm,
97 points
98 );
99
100 pointNormals[pI] += weight * areaNorm;
101 }
102
103 pointNormals[pI].normalise();
105
106 return tpointNormals;
107}
108
109
112(
113 const triSurface& surf,
114 const vectorField& pointNormals
115)
116{
117 const pointField& points = surf.points();
118 const Map<label>& meshPointMap = surf.meshPointMap();
119
120 auto tpointTriads = tmp<triadField>::New(points.size());
121 auto& pointTriads = tpointTriads.ref();
122
123 forAll(points, pI)
124 {
125 const point& pt = points[pI];
126 const vector& normal = pointNormals[meshPointMap[pI]];
127
128 if (mag(normal) < SMALL)
129 {
130 pointTriads[meshPointMap[pI]] = triad::unset;
131 continue;
132 }
133
134 plane p(pt, normal);
135
136 // Pick arbitrary point in plane
137 vector dir1 = normalised(pt - p.somePointInPlane(1e-3));
138 vector dir2 = normalised(dir1 ^ normal);
139
140 pointTriads[meshPointMap[pI]] = triad(dir1, dir2, normal);
142
143 return tpointTriads;
144}
145
146
149(
150 const triSurface& surf,
151 const vectorField& pointNormals,
152 const triadField& pointTriads
153)
154{
155 Info<< "Calculating face curvature" << endl;
156
157 const pointField& points = surf.points();
158 const labelList& meshPoints = surf.meshPoints();
159 const Map<label>& meshPointMap = surf.meshPointMap();
160
161 List<symmTensor2D> pointFundamentalTensors
162 (
163 points.size(),
165 );
166
167 scalarList accumulatedWeights(points.size(), Zero);
168
169 forAll(surf, fI)
170 {
171 const triFace& f = surf[fI];
172 const edgeList fEdges = f.edges();
173
174 // Calculate the edge vectors and the normal differences
175 vectorField edgeVectors(f.size(), Zero);
176 vectorField normalDifferences(f.size(), Zero);
177
178 forAll(fEdges, feI)
179 {
180 const edge& e = fEdges[feI];
181
182 edgeVectors[feI] = e.vec(points);
183 normalDifferences[feI] =
184 pointNormals[meshPointMap[e[0]]]
185 - pointNormals[meshPointMap[e[1]]];
186 }
187
188 // Set up a local coordinate system for the face
189 const vector& e0 = edgeVectors[0];
190 const vector eN = f.areaNormal(points);
191 const vector e1 = (e0 ^ eN);
192
193 if (magSqr(eN) < ROOTVSMALL)
194 {
195 continue;
196 }
197
198 triad faceCoordSys(e0, e1, eN);
199 faceCoordSys.normalize();
200
201 // Construct the matrix to solve
203 scalarDiagonalMatrix Z(3, 0);
204
205 // Least Squares
206 for (label i = 0; i < 3; ++i)
207 {
208 scalar x = edgeVectors[i] & faceCoordSys[0];
209 scalar y = edgeVectors[i] & faceCoordSys[1];
210
211 T(0, 0) += sqr(x);
212 T(1, 0) += x*y;
213 T(1, 1) += sqr(x) + sqr(y);
214 T(2, 1) += x*y;
215 T(2, 2) += sqr(y);
216
217 scalar dndx = normalDifferences[i] & faceCoordSys[0];
218 scalar dndy = normalDifferences[i] & faceCoordSys[1];
219
220 Z[0] += dndx*x;
221 Z[1] += dndx*y + dndy*x;
222 Z[2] += dndy*y;
223 }
224
225 // Perform Cholesky decomposition and back substitution.
226 // Decomposed matrix is in T and solution is in Z.
227 LUsolve(T, Z);
228 symmTensor2D secondFundamentalTensor(Z[0], Z[1], Z[2]);
229
230 // Loop over the face points adding the contribution of the face
231 // curvature to the points.
232 forAll(f, fpI)
233 {
234 const label patchPointIndex = meshPointMap[f[fpI]];
235
236 const triad& ptCoordSys = pointTriads[patchPointIndex];
237
238 if (!ptCoordSys.set())
239 {
240 continue;
241 }
242
243 // Rotate faceCoordSys to ptCoordSys
244 tensor rotTensor = rotationTensor(ptCoordSys[2], faceCoordSys[2]);
245 triad rotatedFaceCoordSys = rotTensor & tensor(faceCoordSys);
246
247 // Project the face curvature onto the point plane
248
249 vector2D cmp1
250 (
251 ptCoordSys[0] & rotatedFaceCoordSys[0],
252 ptCoordSys[0] & rotatedFaceCoordSys[1]
253 );
254 vector2D cmp2
255 (
256 ptCoordSys[1] & rotatedFaceCoordSys[0],
257 ptCoordSys[1] & rotatedFaceCoordSys[1]
258 );
259
260 tensor2D projTensor
261 (
262 cmp1,
263 cmp2
264 );
265
266 symmTensor2D projectedFundamentalTensor
267 (
268 projTensor.x() & (secondFundamentalTensor & projTensor.x()),
269 projTensor.x() & (secondFundamentalTensor & projTensor.y()),
270 projTensor.y() & (secondFundamentalTensor & projTensor.y())
271 );
272
273 // Calculate weight
274 // TODO: Voronoi area weighting
276 (
277 f,
278 meshPoints[patchPointIndex],
279 f.areaNormal(points),
280 points
281 );
282
283 // Sum contribution of face to this point
284 pointFundamentalTensors[patchPointIndex] +=
285 weight*projectedFundamentalTensor;
286
287 accumulatedWeights[patchPointIndex] += weight;
288 }
289
290 if (false)
291 {
292 Info<< "Points = " << points[f[0]] << " "
293 << points[f[1]] << " "
294 << points[f[2]] << endl;
295 Info<< "edgeVecs = " << edgeVectors[0] << " "
296 << edgeVectors[1] << " "
297 << edgeVectors[2] << endl;
298 Info<< "normDiff = " << normalDifferences[0] << " "
299 << normalDifferences[1] << " "
300 << normalDifferences[2] << endl;
301 Info<< "faceCoordSys = " << faceCoordSys << endl;
302 Info<< "T = " << T << endl;
303 Info<< "Z = " << Z << endl;
304 }
305 }
306
307 auto tcurvatureAtPoints = tmp<scalarField>::New(points.size(), Zero);
308 scalarField& curvatureAtPoints = tcurvatureAtPoints.ref();
309
310 forAll(curvatureAtPoints, pI)
311 {
312 pointFundamentalTensors[pI] /= (accumulatedWeights[pI] + SMALL);
313
314 vector2D principalCurvatures(eigenValues(pointFundamentalTensors[pI]));
315
316 //scalar curvature =
317 // (principalCurvatures[0] + principalCurvatures[1])/2;
318 scalar curvature = max
319 (
320 mag(principalCurvatures[0]),
321 mag(principalCurvatures[1])
322 );
323 //scalar curvature = principalCurvatures[0]*principalCurvatures[1];
324
325 curvatureAtPoints[meshPoints[pI]] = curvature;
327
328 return tcurvatureAtPoints;
329}
330
331
334(
335 const triSurface& surf
336)
337{
338 tmp<vectorField> norms = triSurfaceTools::vertexNormals(surf);
339 tmp<triadField> triads = triSurfaceTools::vertexTriads(surf, norms());
340
341 tmp<scalarField> curv = curvatures(surf, norms(), triads());
342 norms.clear();
343 triads.clear();
344
345 return curv;
346}
347
348
351(
352 const Time& runTime,
353 const word& basename,
354 const triSurface& surf
355)
356{
357 Info<< "Extracting curvature of surface at the points." << endl;
358
360 scalarField& curv = tcurv.ref();
361
363 (
365 (
366 basename + ".curvature",
367 runTime.constant(),
368 "triSurface",
369 runTime,
372 ),
373 surf,
374 dimLength,
376 );
377
378 outputField.swap(curv);
379 outputField.write();
380 outputField.swap(curv);
381
382 return tcurv;
383}
384
385
386// ************************************************************************* //
scalar y
void swap(List< T > &other)
Swap with plain List content. Implies shrink_to_fit().
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
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.
Vector2D< Cmpt > y() const
Extract vector for row 1.
Definition Tensor2DI.H:97
Vector2D< Cmpt > x() const
Extract vector for row 0.
Definition Tensor2DI.H:91
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form zero
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
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:74
Representation of a 3D Cartesian coordinate system as a Vector of row vectors.
Definition triad.H:63
static const triad unset
Definition triad.H:107
bool set(const direction d) const
Is the vector in the direction d set.
Definition triadI.H:63
void normalize()
Same as normalise.
Definition triad.H:215
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
const volScalarField & T
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition vector2D.H:56
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedVector eigenValues(const dimensionedSymmTensor &dt)
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition transform.H:47
SymmTensor2D< scalar > symmTensor2D
SymmTensor2D of scalars, i.e. SymmTensor2D<scalar>.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Tensor2D< scalar > tensor2D
Tensor2D of scalars, i.e. Tensor2D<scalar>.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
SymmetricSquareMatrix< scalar > scalarSymmetricSquareMatrix
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
DiagonalMatrix< scalar > scalarDiagonalMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< triad > triadField
Specialisation of Field<T> for triad.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting returning the LU form of the matrix and the sol...
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
DimensionedField< scalar, triSurfacePointGeoMesh > triSurfacePointScalarField
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];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } 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};constexpr 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< DynamicList< face > > pFaces[nBCs]
face triFace(3)
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
3D tensor transformation operations.
Fields for triSurface.