Loading...
Searching...
No Matches
geometricElementTransformPointSmoother.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) 2024-2025 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
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace pointSmoothers
37{
40 (
44 );
45}
46}
47
48
49// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50
51Foam::scalar
52Foam::pointSmoothers::geometricElementTransformPointSmoother::cellQuality
53(
54 const polyMesh& mesh,
55 const pointField& currentPoints,
56 const cellPointConnectivity& connectivity,
57 const label celli
58)
59{
60 const cell& cFaces = mesh.cells()[celli];
61 const labelList cPoints(cFaces.labels(mesh.faces()));
62
63 // Calculate a transformed point for each cell point
64
65 scalar cellQ = 0.0;
66 label nTets = 0;
67
68 forAll(cPoints, cPointi)
69 {
70 const label pointi(cPoints[cPointi]);
71 const point& pt = currentPoints[pointi];
72 const labelList& pPoints
73 (
74 connectivity.cellPointPoints()[celli][cPointi]
75 );
76 if (pPoints.size() != 3)
77 {
78 WarningInFunction<< "Cell:" << celli
79 << " point:" << pointi
80 << " connected points:" << pPoints << endl;
81 }
82 else
83 {
84 const Tensor<scalar> mA
85 (
86 currentPoints[pPoints[0]] - pt,
87 currentPoints[pPoints[1]] - pt,
88 currentPoints[pPoints[2]] - pt
89 );
90 const scalar sigma(det(mA));
91
92 if (sigma < ROOTVSMALL)
93 {
94 return 0;
95 }
96
97 // 3 * pow(sigma, 2.0/3.0)/magSqr(mA)
98 const scalar tetQ =
99 scalar(3) * Foam::cbrt(Foam::sqr(sigma)) / magSqr(mA);
100 cellQ += tetQ;
101 nTets++;
102 }
103 }
105 return cellQ/nTets;
106}
107
108
109// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110
113(
114 const polyMesh& mesh,
115 const dictionary& dict
116)
117:
119 transformationParameter_
120 (
121 dict.get<scalar>("transformationParameter")
122 )
123{}
124
125
126// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
127
129(
130 const labelList& facesToMove,
131 const pointField& oldPoints,
132 const pointField& currentPoints,
133 const pointField& faceCentres,
134 const vectorField& faceAreas,
135 const pointField& cellCentres,
136 const scalarField& cellVolumes,
137 vectorField& pointDisplacement
138) const
139{
140 // Lookup or generate the cell-point connectivity/
141 const cellPointConnectivity& connectivity =
143 (
144 mesh()
145 );
146
147 // Number of points used in each average
148 labelField counts(mesh().nPoints(), -1);
149
150 // Reset the displacements which are about to be calculated
151 reset(facesToMove, counts, pointDisplacement);
152
153 // Identify the cells which are to be moved
154 labelHashSet cellsToMove(facesToMove.size()*2/3);
155 forAll(facesToMove, faceToMoveI)
156 {
157 const label faceI(facesToMove[faceToMoveI]);
158
159 cellsToMove.insert(mesh().faceOwner()[faceI]);
160
161 if (mesh().isInternalFace(faceI))
162 {
163 cellsToMove.insert(mesh().faceNeighbour()[faceI]);
164 }
165 }
166
167 // Transformed point field
168 pointField transformedPoints(currentPoints);
169
170 // Calculate the internal transformations
171 for (const label cellI : cellsToMove)
172 {
173 const cell& cFaces
174 (
175 mesh().cells()[cellI]
176 );
177 const labelList cPoints
178 (
179 cFaces.labels(mesh().faces())
180 );
181 const edgeList cEdges
182 (
183 cFaces.edges(mesh().faces())
184 );
185
186 // Calculate a transformed point for each cell point
187 forAll(cPoints, cPointI)
188 {
189 const label pointI(cPoints[cPointI]);
190
191 if (counts[pointI] == -1) continue;
192
193 const labelList& pPoints
194 (
195 connectivity.cellPointPoints()[cellI][cPointI]
196 );
197 const labelList& pFaces
198 (
199 connectivity.cellPointFaces()[cellI][cPointI]
200 );
201 const label nPPoints(pPoints.size());
202
203 // Initial guess of the dual face centre
204 vector dualAverage(vector::zero);
205 forAll(pPoints, pPointI)
206 {
207 dualAverage +=
208 currentPoints[pPoints[pPointI]]
209 + faceCentres[pFaces[pPointI]];
210 }
211 dualAverage /= 2*nPPoints;
212
213 // Calculate the dual face centre and normal
214 vector dualNormal(vector::zero);
215 forAll(pPoints, pPointI)
216 {
217 const label nextPPointI((pPointI + 1) % nPPoints);
218
219 point edgeCentre
220 (
221 0.5
222 *(
223 currentPoints[pPoints[pPointI]]
224 + currentPoints[pointI]
225 )
226 );
227 point nextFaceCentre
228 (
229 faceCentres[pFaces[nextPPointI]]
230 );
231 point nextEdgeCentre
232 (
233 0.5
234 *(
235 currentPoints[pPoints[nextPPointI]]
236 + currentPoints[pointI]
237 )
238 );
239
240 dualNormal +=
241 (nextFaceCentre - edgeCentre)
242 ^(edgeCentre - dualAverage);
243 dualNormal +=
244 (nextEdgeCentre - nextFaceCentre)
245 ^(nextFaceCentre - dualAverage);
246 }
247 vector dualNormalHat(dualNormal/max(mag(dualNormal), ROOTVSMALL));
248
249 scalar sumA(0);
250 vector sumAc(vector::zero);
251 forAll(pPoints, pPointI)
252 {
253 const label nextPPointI((pPointI + 1) % nPPoints);
254
255 point edgeCentre
256 (
257 0.5
258 *(
259 currentPoints[pPoints[pPointI]]
260 + currentPoints[pointI]
261 )
262 );
263 point nextFaceCentre
264 (
265 faceCentres[pFaces[nextPPointI]]
266 );
267 point nextEdgeCentre
268 (
269 0.5
270 *(
271 currentPoints[pPoints[nextPPointI]]
272 + currentPoints[pointI]
273 )
274 );
275
276 vector c1 = edgeCentre + nextFaceCentre + dualAverage;
277 vector c2 = nextFaceCentre + nextEdgeCentre + dualAverage;
278
279 vector n1 =
280 (nextFaceCentre - edgeCentre)
281 ^(edgeCentre - dualAverage);
282 vector n2 =
283 (nextEdgeCentre - nextFaceCentre)
284 ^(nextFaceCentre - dualAverage);
285
286 scalar a1 = n1 & dualNormalHat;
287 scalar a2 = n2 & dualNormalHat;
288
289 sumA += a1 + a2;
290 sumAc += a1*c1 + a2*c2;
291 }
292
293 const vector dualCentre(sumAc/max(sumA, ROOTVSMALL)/3);
294
295 // Calculate the transformed point
296 transformedPoints[pointI] =
297 dualCentre
298 + transformationParameter_
299 *dualNormal/sqrt(max(mag(dualNormal), ROOTVSMALL));
300 }
301
302 // Length scale
303 scalar lengthScale(0), transformedLengthScale(0);
304 forAll(cEdges, cEdgeI)
305 {
306 lengthScale +=
307 cEdges[cEdgeI].mag(currentPoints);
308 transformedLengthScale +=
309 cEdges[cEdgeI].mag(transformedPoints);
310 }
311 lengthScale /= cEdges.size();
312 transformedLengthScale /= cEdges.size();
313
314 const scalar lengthScaleRatio =
315 (
316 (transformedLengthScale > SMALL)
317 ? lengthScale/transformedLengthScale
318 : scalar(0)
319 );
320
321 // Add the displacement to the average
322 forAll(cPoints, cPointI)
323 {
324 const label pointI(cPoints[cPointI]);
325
326 if (counts[pointI] == -1) continue;
327
328 const vector newPoint
329 (
330 cellCentres[cellI]
331 + lengthScaleRatio
332 *(
333 transformedPoints[pointI]
334 - cellCentres[cellI]
335 )
336 );
337
338 ++ counts[pointI];
339
340 pointDisplacement[pointI] += newPoint - oldPoints[pointI];
341 }
342 }
343
344 // Reset all the boundary faces
345 reset(facesToMove, counts, pointDisplacement, false);
346
347 // Calculate the boundary transformations
348 forAll(facesToMove, faceToMoveI)
349 {
350 const label faceI(facesToMove[faceToMoveI]);
351
352 if (!isInternalOrProcessorFace(faceI))
353 {
354 const labelList& fPoints(mesh().faces()[faceI]);
355
356 // Face normal
357 vector faceNormalHat(faceAreas[faceI]);
358 faceNormalHat /= max(SMALL, mag(faceNormalHat));
359
360 // Calculate a transformed point for each face point
361 forAll(fPoints, fPointI)
362 {
363 const label pointI(fPoints[fPointI]);
364
365 const label fPointIPrev
366 (
367 (fPointI - 1 + fPoints.size()) % fPoints.size()
368 );
369 const label fPointINext
370 (
371 (fPointI + 1) % fPoints.size()
372 );
373
374 const vector dualCentre
375 (
376 currentPoints[pointI]/2
377 + (
378 currentPoints[fPoints[fPointINext]]
379 + currentPoints[fPoints[fPointIPrev]]
380 )/4
381 );
382 const vector dualNormal
383 (
384 (
385 currentPoints[fPoints[fPointINext]]
386 - currentPoints[fPoints[fPointIPrev]]
387 )/2
388 ^faceNormalHat
389 );
390
391 transformedPoints[pointI] =
392 dualCentre
393 + transformationParameter_
394 *dualNormal/sqrt(max(mag(dualNormal), ROOTVSMALL));
395 }
396
397 // Length scale
398 scalar lengthScale(0), transformedLengthScale(0);
399 forAll(fPoints, fPointI)
400 {
401 const label fPointINext((fPointI + 1)%fPoints.size());
402
403 const label pointI(fPoints[fPointI]);
404 const label pointINext(fPoints[fPointINext]);
405
406 lengthScale += mag
407 (
408 currentPoints[pointINext] - currentPoints[pointI]
409 );
410 transformedLengthScale += mag
411 (
412 transformedPoints[pointINext] - transformedPoints[pointI]
413 );
414 }
415 lengthScale /= fPoints.size();
416 transformedLengthScale /= fPoints.size();
417
418 const scalar lengthScaleRatio
419 (
420 (transformedLengthScale > SMALL)
421 ? lengthScale/transformedLengthScale
422 : scalar(0)
423 );
424
425 // Add the displacement to the average
426 forAll(fPoints, fPointI)
427 {
428 const label pointI(fPoints[fPointI]);
429
430 const vector newPoint
431 (
432 faceCentres[faceI]
433 + lengthScaleRatio
434 *(
435 transformedPoints[pointI]
436 - faceCentres[faceI]
437 )
438 );
439
440 ++ counts[pointI];
441
442 pointDisplacement[pointI] += newPoint - oldPoints[pointI];
443 }
444 }
445 }
447 // Average
448 average(facesToMove, counts, pointDisplacement);
449}
450
451
453Foam::pointSmoothers::geometricElementTransformPointSmoother::cellQuality
454(
455 const polyMesh& mesh,
456 const pointField& currentPoints
457)
458{
459 const cellPointConnectivity& connectivity =
461 (
462 mesh
463 );
464
466 scalarField& fld = tfld.ref();
467
468 forAll(fld, celli)
469 {
470 fld[celli] = cellQuality(mesh, currentPoints, connectivity, celli);
471 }
472
473 return tfld;
474}
475
476
477// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
static FOAM_NO_DANGLING_REFERENCE const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form zero
This class provides ordered connectivity for each point of each cell. Lists are available of the poin...
const labelListListList & cellPointFaces() const
Access the point-face connections.
const labelListListList & cellPointPoints() const
Access the point-point connections.
Class calculates cell quality measures.
Definition cellQuality.H:48
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
labelList labels(const faceUList &meshFaces) const
Return unordered list of cell vertices given the list of faces.
Definition cell.C:33
edgeList edges(const faceUList &meshFaces) const
Return cell edges.
Definition cell.C:105
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for point smoothing methods. Handles parallel communication via reset and average...
bool isInternalOrProcessorFace(const label faceI) const
Test if the given face is internal or on a processor boundary.
const polyMesh & mesh() const noexcept
Access the mesh.
Geometric Element Transformation Method (GETMe) point smoother. Points are moved in a direction norma...
virtual void calculate(const labelList &facesToMove, const pointField &oldPoints, const pointField &currentPoints, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes, vectorField &pointDisplacement) const
Calculate the point displacements.
geometricElementTransformPointSmoother(const polyMesh &mesh, const dictionary &dict)
Construct from a dictionary and a mesh.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A class for managing temporary objects.
Definition tmp.H:75
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 Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
limits reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL))
dynamicFvMesh & mesh
label nPoints
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
dimensionedScalar det(const dimensionedSphericalTensor &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
vector point
Point is a vector.
Definition point.H:37
dimensionedScalar cbrt(const dimensionedScalar &ds)
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
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];} 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]
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299