polyMeshTetDecomposition.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) 2011-2019 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 // Note: the use of this tolerance is ad-hoc, there may be extreme
34 // cases where the resulting tetrahedra still have particle tracking
35 // problems, or tets with lower quality may track OK.
36 const Foam::scalar Foam::polyMeshTetDecomposition::minTetQuality = sqr(SMALL);
37 
38 
39 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
40 
42 (
43  const polyMesh& mesh,
44  const point& cC,
45  const label fI,
46  const bool isOwner,
47  const label faceBasePtI
48 )
49 {
50  // Does fan decomposition of face (starting at faceBasePti) and determines
51  // min quality over all resulting tets.
52 
53  const pointField& pPts = mesh.points();
54  const face& f = mesh.faces()[fI];
55  const point& tetBasePt = pPts[f[faceBasePtI]];
56 
57  scalar thisBaseMinTetQuality = VGREAT;
58 
59  for (label tetPtI = 1; tetPtI < f.size() - 1; tetPtI++)
60  {
61  label facePtI = (tetPtI + faceBasePtI) % f.size();
62  label otherFacePtI = f.fcIndex(facePtI);
63 
64  label ptAI = -1;
65  label ptBI = -1;
66 
67  if (isOwner)
68  {
69  ptAI = f[facePtI];
70  ptBI = f[otherFacePtI];
71  }
72  else
73  {
74  ptAI = f[otherFacePtI];
75  ptBI = f[facePtI];
76  }
77 
78  const point& pA = pPts[ptAI];
79  const point& pB = pPts[ptBI];
80 
81  tetPointRef tet(cC, tetBasePt, pA, pB);
82 
83  scalar tetQuality = tet.quality();
84 
85  if (tetQuality < thisBaseMinTetQuality)
86  {
87  thisBaseMinTetQuality = tetQuality;
88  }
89  }
90  return thisBaseMinTetQuality;
91 }
92 
93 
95 (
96  const polyMesh& mesh,
97  label fI,
98  const point& nCc,
99  scalar tol,
100  bool report
101 )
102 {
103  const faceList& pFaces = mesh.faces();
104  const vectorField& pC = mesh.cellCentres();
105  const labelList& pOwner = mesh.faceOwner();
106 
107  const face& f = pFaces[fI];
108 
109  label oCI = pOwner[fI];
110 
111  const point& oCc = pC[oCI];
112 
113  forAll(f, faceBasePtI)
114  {
115  scalar minQ = minQuality(mesh, oCc, fI, true, faceBasePtI);
116  minQ = min(minQ, minQuality(mesh, nCc, fI, false, faceBasePtI));
117 
118  if (minQ > tol)
119  {
120  return faceBasePtI;
121  }
122  }
123 
124  // If a base point hasn't triggered a return by now, then there is
125  // none that can produce a good decomposition
126  return -1;
127 }
128 
129 
131 (
132  const polyMesh& mesh,
133  label fI,
134  scalar tol,
135  bool report
136 )
137 {
138  return findSharedBasePoint
139  (
140  mesh,
141  fI,
142  mesh.cellCentres()[mesh.faceNeighbour()[fI]],
143  tol,
144  report
145  );
146 }
147 
148 
150 (
151  const polyMesh& mesh,
152  label fI,
153  scalar tol,
154  bool report
155 )
156 {
157  const faceList& pFaces = mesh.faces();
158  const vectorField& pC = mesh.cellCentres();
159  const labelList& pOwner = mesh.faceOwner();
160 
161  const face& f = pFaces[fI];
162 
163  label cI = pOwner[fI];
164 
165  bool own = (pOwner[fI] == cI);
166 
167  const point& cC = pC[cI];
168 
169  forAll(f, faceBasePtI)
170  {
171  scalar quality = minQuality(mesh, cC, fI, own, faceBasePtI);
172 
173  if (quality > tol)
174  {
175  return faceBasePtI;
176  }
177  }
178 
179  // If a base point hasn't triggered a return by now, then there is
180  // none that can produce a good decomposition
181  return -1;
182 }
183 
184 
186 (
187  const polyMesh& mesh,
188  scalar tol,
189  bool report
190 )
191 {
192  const labelList& pOwner = mesh.faceOwner();
193  const vectorField& pC = mesh.cellCentres();
194 
195  // Find a suitable base point for each face, considering both
196  // cells for interface faces or those on coupled patches
197 
198  labelList tetBasePtIs(mesh.nFaces(), -1);
199 
200  label nInternalFaces = mesh.nInternalFaces();
201 
202  for (label fI = 0; fI < nInternalFaces; ++fI)
203  {
204  tetBasePtIs[fI] = findSharedBasePoint(mesh, fI, tol, report);
205  }
206 
207  pointField neighbourCellCentres(mesh.nFaces() - nInternalFaces);
208 
209  for (label facei = nInternalFaces; facei < mesh.nFaces(); ++facei)
210  {
211  neighbourCellCentres[facei - nInternalFaces] = pC[pOwner[facei]];
212  }
213 
214  syncTools::swapBoundaryFacePositions(mesh, neighbourCellCentres);
215 
216  const polyBoundaryMesh& patches = mesh.boundaryMesh();
217 
218  SubList<label> boundaryFaceTetBasePtIs
219  (
220  tetBasePtIs,
221  mesh.nFaces() - nInternalFaces,
222  nInternalFaces
223  );
224 
225  for
226  (
227  label fI = nInternalFaces, bFI = 0;
228  fI < mesh.nFaces();
229  fI++, bFI++
230  )
231  {
232  label patchi = mesh.boundaryMesh().patchID()[bFI];
233 
234  if (patches[patchi].coupled())
235  {
236  const coupledPolyPatch& pp =
237  refCast<const coupledPolyPatch>(patches[patchi]);
238 
239  if (pp.owner())
240  {
241  boundaryFaceTetBasePtIs[bFI] = findSharedBasePoint
242  (
243  mesh,
244  fI,
245  neighbourCellCentres[bFI],
246  tol,
247  report
248  );
249  }
250  else
251  {
252  // Assign -2, to distinguish from a failed base point
253  // find, which returns -1.
254  boundaryFaceTetBasePtIs[bFI] = -2;
255  }
256  }
257  else
258  {
259  boundaryFaceTetBasePtIs[bFI] = findBasePoint
260  (
261  mesh,
262  fI,
263  tol,
264  report
265  );
266  }
267  }
268 
269  // maxEqOp will replace the -2 values on the neighbour patches
270  // with the result from the owner base point find.
271 
272  syncTools::syncBoundaryFaceList
273  (
274  mesh,
275  boundaryFaceTetBasePtIs,
277  );
278 
279  for
280  (
281  label fI = nInternalFaces, bFI = 0;
282  fI < mesh.nFaces();
283  fI++, bFI++
284  )
285  {
286  label& bFTetBasePtI = boundaryFaceTetBasePtIs[bFI];
287 
288  if (bFTetBasePtI == -2)
289  {
291  << "Coupled face base point exchange failure for face "
292  << fI << " at " << mesh.faceCentres()[fI]
293  << abort(FatalError);
294  }
295 
296  if (bFTetBasePtI < 1)
297  {
298  // If the base point is -1, it should be left as such to
299  // indicate a problem, if it is 0, then no action is required.
300 
301  continue;
302  }
303 
304  label patchi = mesh.boundaryMesh().patchID()[bFI];
305 
306  if (patches[patchi].coupled())
307  {
308  const coupledPolyPatch& pp =
309  refCast<const coupledPolyPatch>(patches[patchi]);
310 
311  // Calculated base points on coupled faces are those of
312  // the owner patch face. They need to be reindexed to for
313  // the non-owner face, which has the opposite order.
314 
315  // So, for fPtI_o != 0, fPtI_n = f.size() - fPtI_o
316 
317  // i.e.:
318 
319  // owner coupledPolyPatch face
320  // face (a b c d e f)
321  // fPtI 0 1 2 3 4 5
322  // +
323  // #
324 
325  // neighbour coupledPolyPatch face
326  // face (a f e d c b)
327  // fPtI 0 1 2 3 4 5
328  // +
329  // #
330  // +: 6 - 1 = 5
331  // #: 6 - 2 = 4
332 
333  if (!pp.owner())
334  {
335  bFTetBasePtI = mesh.faces()[fI].size() - bFTetBasePtI;
336  }
337  }
338  }
339 
340  return tetBasePtIs;
341 }
342 
343 
345 (
346  const polyMesh& mesh,
347  scalar tol,
348  const bool report,
349  labelHashSet* setPtr
350 )
351 {
352  const labelList& own = mesh.faceOwner();
353  const labelList& nei = mesh.faceNeighbour();
354  const polyBoundaryMesh& patches = mesh.boundaryMesh();
355 
356  const vectorField& cc = mesh.cellCentres();
357  const vectorField& fc = mesh.faceCentres();
358 
359  // Calculate coupled cell centre
360  pointField neiCc(mesh.nBoundaryFaces());
361 
362  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
363  {
364  neiCc[facei - mesh.nInternalFaces()] = cc[own[facei]];
365  }
366 
367  syncTools::swapBoundaryFacePositions(mesh, neiCc);
368 
369  const faceList& fcs = mesh.faces();
370 
371  const pointField& p = mesh.points();
372 
373  label nErrorTets = 0;
374 
375  forAll(fcs, facei)
376  {
377  const face& f = fcs[facei];
378 
379  forAll(f, fPtI)
380  {
381  scalar tetQual = tetPointRef
382  (
383  p[f[fPtI]],
384  p[f.nextLabel(fPtI)],
385  fc[facei],
386  cc[own[facei]]
387  ).quality();
388 
389  if (tetQual > -tol)
390  {
391  if (setPtr)
392  {
393  setPtr->insert(facei);
394  }
395 
396  nErrorTets++;
397  break; // no need to check other tets
398  }
399  }
400 
401  if (mesh.isInternalFace(facei))
402  {
403  // Create the neighbour tet - it will have positive volume
404  const face& f = fcs[facei];
405 
406  forAll(f, fPtI)
407  {
408  scalar tetQual = tetPointRef
409  (
410  p[f[fPtI]],
411  p[f.nextLabel(fPtI)],
412  fc[facei],
413  cc[nei[facei]]
414  ).quality();
415 
416  if (tetQual < tol)
417  {
418  if (setPtr)
419  {
420  setPtr->insert(facei);
421  }
422 
423  nErrorTets++;
424  break;
425  }
426  }
427 
428  if (findSharedBasePoint(mesh, facei, tol, report) == -1)
429  {
430  if (setPtr)
431  {
432  setPtr->insert(facei);
433  }
434 
435  nErrorTets++;
436  }
437  }
438  else
439  {
440  label patchi = patches.patchID()[facei - mesh.nInternalFaces()];
441 
442  if (patches[patchi].coupled())
443  {
444  if
445  (
446  findSharedBasePoint
447  (
448  mesh,
449  facei,
450  neiCc[facei - mesh.nInternalFaces()],
451  tol,
452  report
453  ) == -1
454  )
455  {
456  if (setPtr)
457  {
458  setPtr->insert(facei);
459  }
460 
461  nErrorTets++;
462  }
463  }
464  else
465  {
466  if (findBasePoint(mesh, facei, tol, report) == -1)
467  {
468  if (setPtr)
469  {
470  setPtr->insert(facei);
471  }
472 
473  nErrorTets++;
474  }
475  }
476  }
477  }
478 
479  reduce(nErrorTets, sumOp<label>());
480 
481  if (nErrorTets > 0)
482  {
483  if (report)
484  {
485  Info<< " ***Error in face tets: "
486  << nErrorTets << " faces with low quality or negative volume "
487  << "decomposition tets." << endl;
488  }
489 
490  return true;
491  }
492 
493  if (report)
494  {
495  Info<< " Face tets OK." << endl;
496  }
497 
498  return false;
499 }
500 
501 
503 (
504  const polyMesh& mesh,
505  label facei,
506  label celli
507 )
508 {
509  const faceList& pFaces = mesh.faces();
510 
511  const face& f = pFaces[facei];
512 
513  label nTets = f.size() - 2;
514 
515  List<tetIndices> faceTets(nTets);
516 
517  for (label tetPti = 1; tetPti < f.size() - 1; ++tetPti)
518  {
519  faceTets[tetPti - 1] = tetIndices(celli, facei, tetPti);
520  }
521 
522  return faceTets;
523 }
524 
525 
527 (
528  const polyMesh& mesh,
529  label celli
530 )
531 {
532  const faceList& pFaces = mesh.faces();
533  const cellList& pCells = mesh.cells();
534 
535  const cell& thisCell = pCells[celli];
536 
537  label nTets = 0;
538 
539  for (const label facei : thisCell)
540  {
541  nTets += pFaces[facei].size() - 2;
542  }
543 
544  DynamicList<tetIndices> cellTets(nTets);
545 
546  for (const label facei : thisCell)
547  {
548  cellTets.append(faceTetIndices(mesh, facei, celli));
549  }
550 
551  return cellTets;
552 }
553 
554 
556 (
557  const polyMesh& mesh,
558  label celli,
559  const point& pt
560 )
561 {
562  const faceList& pFaces = mesh.faces();
563  const cellList& pCells = mesh.cells();
564 
565  const cell& thisCell = pCells[celli];
566 
567  tetIndices tetContainingPt;
568 
569 
570  for (const label facei : thisCell)
571  {
572  const face& f = pFaces[facei];
573 
574  for (label tetPti = 1; tetPti < f.size() - 1; ++tetPti)
575  {
576  // Get tetIndices of face triangle
577  tetIndices faceTetIs(celli, facei, tetPti);
578 
579  // Check if inside
580  if (faceTetIs.tet(mesh).inside(pt))
581  {
582  tetContainingPt = faceTetIs;
583  break;
584  }
585  }
586 
587  if (tetContainingPt.cell() != -1)
588  {
589  break;
590  }
591  }
592 
593  return tetContainingPt;
594 }
595 
596 
597 // ************************************************************************* //
Foam::polyMeshTetDecomposition::checkFaceTets
static bool checkFaceTets(const polyMesh &mesh, scalar tol=minTetQuality, const bool report=false, labelHashSet *setPtr=nullptr)
Check face-decomposition tet volume.
Definition: polyMeshTetDecomposition.C:345
Foam::tetrahedron::inside
bool inside(const point &pt) const
Return true if point is inside tetrahedron.
Definition: tetrahedronI.H:415
Foam::polyMeshTetDecomposition::minTetQuality
static const scalar minTetQuality
Minimum tetrahedron quality.
Definition: polyMeshTetDecomposition.H:66
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::tetIndices::tet
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:154
Foam::tetIndices::cell
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::polyMeshTetDecomposition::findFaceBasePts
static labelList findFaceBasePts(const polyMesh &mesh, scalar tol=minTetQuality, bool report=false)
Find a suitable base point for each face for decomposition.
Definition: polyMeshTetDecomposition.C:186
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::coupledPolyPatch
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Definition: coupledPolyPatch.H:53
Foam::coupledPolyPatch::owner
virtual bool owner() const =0
Does this side own the patch ?
Foam::HashSet< label, Hash< label > >
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
polyMeshTetDecomposition.H
Foam::sumOp
Definition: ops.H:213
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:47
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::polyMeshTetDecomposition::findSharedBasePoint
static label findSharedBasePoint(const polyMesh &mesh, label fI, const point &nCc, scalar tol, bool report=false)
Find the first suitable base point to use for a minimum.
Definition: polyMeshTetDecomposition.C:95
Foam::polyMeshTetDecomposition::findBasePoint
static label findBasePoint(const polyMesh &mesh, label fI, scalar tol, bool report=false)
Find the base point to use for a minimum triangle.
Definition: polyMeshTetDecomposition.C:150
pFaces
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]=cellShape(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
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:474
Foam::tetrahedron::quality
scalar quality() const
Return quality: Ratio of tetrahedron and circum-sphere.
Definition: tetrahedronI.H:232
Foam::FatalError
error FatalError
reduce
reduce(hasMovingMesh, orOp< bool >())
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::polyMeshTetDecomposition::findTet
static tetIndices findTet(const polyMesh &mesh, label cI, const point &pt)
Find the tet decomposition of the cell containing the given point.
Definition: polyMeshTetDecomposition.C:556
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::maxEqOp
Definition: ops.H:80
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< face >
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:181
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::polyMeshTetDecomposition::cellTetIndices
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
Definition: polyMeshTetDecomposition.C:527
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::polyMeshTetDecomposition::minQuality
static scalar minQuality(const polyMesh &mesh, const point &cC, const label fI, const bool isOwner, const label faceBasePtI)
Given a face and cc and starting index for triangulation determine.
Definition: polyMeshTetDecomposition.C:42
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:65
Foam::polyMeshTetDecomposition::faceTetIndices
static List< tetIndices > faceTetIndices(const polyMesh &mesh, label fI, label cI)
Return the tet decomposition of the given face, with.
Definition: polyMeshTetDecomposition.C:503