polyMeshGeometry.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-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2022 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 
29 #include "polyMeshGeometry.H"
31 #include "pyramid.H"
32 #include "tetrahedron.H"
33 #include "syncTools.H"
34 #include "unitConversion.H"
35 #include "primitiveMeshTools.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(polyMeshGeometry, 0);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::polyMeshGeometry::updateFaceCentresAndAreas
48 (
49  const pointField& p,
50  const labelList& changedFaces
51 )
52 {
53  const faceList& fcs = mesh_.faces();
54 
55  for (const label facei : changedFaces)
56  {
57  const face& f = fcs[facei];
58  const label nPoints = f.size();
59 
60  // If the face is a triangle, do a direct calculation for efficiency
61  // and to avoid round-off error-related problems
62  if (nPoints == 3)
63  {
64  faceCentres_[facei] =
65  triPointRef::centre(p[f[0]], p[f[1]], p[f[2]]);
66  faceAreas_[facei] =
67  triPointRef::areaNormal(p[f[0]], p[f[1]], p[f[2]]);
68  }
69  else
70  {
71  solveVector sumN = Zero;
72  solveScalar sumA = Zero;
73  solveVector sumAc = Zero;
74 
75  solveVector fCentre = p[f[0]];
76  for (label pi = 1; pi < nPoints; ++pi)
77  {
78  fCentre += solveVector(p[f[pi]]);
79  }
80  fCentre /= nPoints;
81 
82  for (label pi = 0; pi < nPoints; ++pi)
83  {
84  const solveVector thisPoint(p[f.thisLabel(pi)]);
85  const solveVector nextPoint(p[f.nextLabel(pi)]);
86 
87  solveVector c = thisPoint + nextPoint + fCentre;
88  solveVector n = (nextPoint - thisPoint)^(fCentre - thisPoint);
89  solveScalar a = mag(n);
90 
91  sumN += n;
92  sumA += a;
93  sumAc += a*c;
94  }
95 
96  if (sumA < ROOTVSMALL)
97  {
98  faceCentres_[facei] = fCentre;
99  faceAreas_[facei] = Zero;
100  }
101  else
102  {
103  faceCentres_[facei] = (1.0/3.0)*sumAc/sumA;
104  faceAreas_[facei] = 0.5*sumN;
105  }
106  }
107  }
108 }
109 
110 
111 void Foam::polyMeshGeometry::updateCellCentresAndVols
112 (
113  const labelList& changedCells,
114  const labelList& changedFaces
115 )
116 {
117  const labelList& own = mesh().faceOwner();
118  const cellList& cells = mesh().cells();
119 
120  // Clear the fields for accumulation
121  UIndirectList<vector>(cellCentres_, changedCells) = Zero;
122  UIndirectList<scalar>(cellVolumes_, changedCells) = Zero;
123 
124 
125  // Re-calculate the changed cell centres and volumes
126  for (const label celli : changedCells)
127  {
128  const labelList& cFaces = cells[celli];
129 
130  // Estimate the cell centre and bounding box using the face centres
131  vector cEst(Zero);
132  boundBox bb;
133 
134  for (const label facei : cFaces)
135  {
136  const point& fc = faceCentres_[facei];
137  cEst += fc;
138  bb.add(fc);
139  }
140  cEst /= cFaces.size();
141 
142 
143  // Sum up the face-pyramid contributions
144  for (const label facei : cFaces)
145  {
146  // Calculate 3* the face-pyramid volume
147  scalar pyr3Vol = faceAreas_[facei] & (faceCentres_[facei] - cEst);
148 
149  if (own[facei] != celli)
150  {
151  pyr3Vol = -pyr3Vol;
152  }
153 
154  // Accumulate face-pyramid volume
155  cellVolumes_[celli] += pyr3Vol;
156 
157  // Calculate the face-pyramid centre
158  const vector pCtr = (3.0/4.0)*faceCentres_[facei] + (1.0/4.0)*cEst;
159 
160  // Accumulate volume-weighted face-pyramid centre
161  cellCentres_[celli] += pyr3Vol*pCtr;
162  }
163 
164  // Average the accumulated quantities
165 
166  if (mag(cellVolumes_[celli]) > VSMALL)
167  {
168  point cc = cellCentres_[celli] / cellVolumes_[celli];
169 
170  // Do additional check for collapsed cells since some volumes
171  // (e.g. 1e-33) do not trigger above but do return completely
172  // wrong cell centre
173  if (bb.contains(cc))
174  {
175  cellCentres_[celli] = cc;
176  }
177  else
178  {
179  cellCentres_[celli] = cEst;
180  }
181  }
182  else
183  {
184  cellCentres_[celli] = cEst;
185  }
187  cellVolumes_[celli] *= (1.0/3.0);
188  }
189 }
190 
191 
193 (
194  const polyMesh& mesh,
195  const labelList& changedFaces
196 )
197 {
198  const labelList& own = mesh.faceOwner();
199  const labelList& nei = mesh.faceNeighbour();
200 
201  labelHashSet affectedCells(2*changedFaces.size());
202 
203  for (const label facei : changedFaces)
204  {
205  affectedCells.insert(own[facei]);
206 
207  if (mesh.isInternalFace(facei))
208  {
209  affectedCells.insert(nei[facei]);
210  }
211  }
212  return affectedCells.toc();
213 }
214 
215 
216 Foam::scalar Foam::polyMeshGeometry::checkNonOrtho
217 (
218  const polyMesh& mesh,
219  const bool report,
220  const scalar severeNonorthogonalityThreshold,
221  const label facei,
222  const vector& s, // face area vector
223  const vector& d, // cc-cc vector
224 
225  label& severeNonOrth,
226  label& errorNonOrth,
227  labelHashSet* setPtr
228 )
229 {
230  scalar dDotS = (d & s)/(mag(d)*mag(s) + VSMALL);
231 
232  if (dDotS < severeNonorthogonalityThreshold)
233  {
234  label nei = -1;
235 
236  if (mesh.isInternalFace(facei))
237  {
238  nei = mesh.faceNeighbour()[facei];
239  }
240 
241  if (dDotS > SMALL)
242  {
243  if (report)
244  {
245  // Severe non-orthogonality but mesh still OK
246  Pout<< "Severe non-orthogonality for face " << facei
247  << " between cells " << mesh.faceOwner()[facei]
248  << " and " << nei
249  << ": Angle = "
250  << radToDeg(::acos(dDotS))
251  << " deg." << endl;
252  }
253 
254  ++severeNonOrth;
255  }
256  else
257  {
258  // Non-orthogonality greater than 90 deg
259  if (report)
260  {
262  << "Severe non-orthogonality detected for face "
263  << facei
264  << " between cells " << mesh.faceOwner()[facei]
265  << " and " << nei
266  << ": Angle = "
267  << radToDeg(::acos(dDotS))
268  << " deg." << endl;
269  }
270 
271  ++errorNonOrth;
272  }
273 
274  if (setPtr)
275  {
276  setPtr->insert(facei);
277  }
278  }
279  return dDotS;
280 }
281 
282 
283 // Create the neighbour pyramid - it will have positive volume
284 bool Foam::polyMeshGeometry::checkFaceTet
285 (
286  const polyMesh& mesh,
287  const bool report,
288  const scalar minTetQuality,
289  const pointField& p,
290  const label facei,
291  const point& fc, // face centre
292  const point& cc, // cell centre
293 
294  labelHashSet* setPtr
295 )
296 {
297  const face& f = mesh.faces()[facei];
298 
299  forAll(f, fp)
300  {
301  scalar tetQual = tetPointRef
302  (
303  p[f[fp]],
304  p[f.nextLabel(fp)],
305  fc,
306  cc
307  ).quality();
308 
309  if (tetQual < minTetQuality)
310  {
311  if (report)
312  {
313  Pout<< "bool polyMeshGeometry::checkFaceTets("
314  << "const bool, const scalar, const pointField&"
315  << ", const pointField&"
316  << ", const labelList&, labelHashSet*) : "
317  << "face " << facei
318  << " has a triangle that points the wrong way." << nl
319  << "Tet quality: " << tetQual
320  << " Face " << facei
321  << endl;
322  }
323  if (setPtr)
324  {
325  setPtr->insert(facei);
326  }
327  return true;
328  }
329  }
331  return false;
332 }
333 
334 
335 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
336 
338 :
339  mesh_(mesh)
340 {
341  correct();
342 }
343 
344 
345 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
346 
348 {
349  faceAreas_ = mesh_.faceAreas();
350  faceCentres_ = mesh_.faceCentres();
351  cellCentres_ = mesh_.cellCentres();
352  cellVolumes_ = mesh_.cellVolumes();
353 }
354 
355 
357 (
358  const pointField& p,
359  const labelList& changedFaces
360 )
361 {
362  // Update face quantities
363  updateFaceCentresAndAreas(p, changedFaces);
364  // Update cell quantities from face quantities
365  updateCellCentresAndVols(affectedCells(mesh_, changedFaces), changedFaces);
366 }
367 
368 
370 (
371  const bool report,
372  const scalar orthWarn,
373  const polyMesh& mesh,
374  const vectorField& cellCentres,
375  const vectorField& faceAreas,
376  const labelList& checkFaces,
377  const List<labelPair>& baffles,
378  labelHashSet* setPtr
379 )
380 {
381  // for all internal and coupled faces check that the d dot S product
382  // is positive
383 
384  const labelList& own = mesh.faceOwner();
385  const labelList& nei = mesh.faceNeighbour();
387 
388  // Severe nonorthogonality threshold
389  const scalar severeNonorthogonalityThreshold = ::cos(degToRad(orthWarn));
390 
391  // Calculate coupled cell centre
392  pointField neiCc(mesh.nBoundaryFaces());
393 
394  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
395  {
396  neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
397  }
398 
400 
401  scalar minDDotS = GREAT;
402 
403  scalar sumDDotS = 0;
404  label nDDotS = 0;
405 
406  label severeNonOrth = 0;
407 
408  label errorNonOrth = 0;
409 
410  for (const label facei : checkFaces)
411  {
412  const point& ownCc = cellCentres[own[facei]];
413 
414  if (mesh.isInternalFace(facei))
415  {
416  scalar dDotS = checkNonOrtho
417  (
418  mesh,
419  report,
420  severeNonorthogonalityThreshold,
421  facei,
422  faceAreas[facei],
423  cellCentres[nei[facei]] - ownCc,
424 
425  severeNonOrth,
426  errorNonOrth,
427  setPtr
428  );
429 
430  if (dDotS < minDDotS)
431  {
432  minDDotS = dDotS;
433  }
434 
435  sumDDotS += dDotS;
436  ++nDDotS;
437  }
438  else
439  {
440  const label patchi = patches.whichPatch(facei);
441 
442  if (patches[patchi].coupled())
443  {
444  scalar dDotS = checkNonOrtho
445  (
446  mesh,
447  report,
448  severeNonorthogonalityThreshold,
449  facei,
450  faceAreas[facei],
451  neiCc[facei-mesh.nInternalFaces()] - ownCc,
452 
453  severeNonOrth,
454  errorNonOrth,
455  setPtr
456  );
457 
458  if (dDotS < minDDotS)
459  {
460  minDDotS = dDotS;
461  }
462 
463  sumDDotS += dDotS;
464  ++nDDotS;
465  }
466  }
467  }
468 
469  for (const labelPair& baffle : baffles)
470  {
471  const label face0 = baffle.first();
472  const label face1 = baffle.second();
473 
474  const point& ownCc = cellCentres[own[face0]];
475 
476  scalar dDotS = checkNonOrtho
477  (
478  mesh,
479  report,
480  severeNonorthogonalityThreshold,
481  face0,
482  faceAreas[face0],
483  cellCentres[own[face1]] - ownCc,
484 
485  severeNonOrth,
486  errorNonOrth,
487  setPtr
488  );
489 
490  if (dDotS < minDDotS)
491  {
492  minDDotS = dDotS;
493  }
494 
495  sumDDotS += dDotS;
496  ++nDDotS;
497  }
498 
499  reduce(minDDotS, minOp<scalar>());
500  reduce(sumDDotS, sumOp<scalar>());
501  reduce(nDDotS, sumOp<label>());
502  reduce(severeNonOrth, sumOp<label>());
503  reduce(errorNonOrth, sumOp<label>());
504 
505  // Only report if there are some internal faces
506  if (nDDotS > 0)
507  {
508  if (report && minDDotS < severeNonorthogonalityThreshold)
509  {
510  Info<< "Number of non-orthogonality errors: " << errorNonOrth
511  << ". Number of severely non-orthogonal faces: "
512  << severeNonOrth << "." << endl;
513  }
514  }
515 
516  if (report)
517  {
518  if (nDDotS > 0)
519  {
520  Info<< "Mesh non-orthogonality Max: "
521  << radToDeg(::acos(minDDotS))
522  << " average: " << radToDeg(::acos(sumDDotS/nDDotS))
523  << endl;
524  }
525  }
526 
527  if (errorNonOrth > 0)
528  {
529  if (report)
530  {
532  << "Error in non-orthogonality detected" << endl;
533  }
534 
535  return true;
536  }
537 
538  if (report)
539  {
540  Info<< "Non-orthogonality check OK.\n" << endl;
541  }
542 
543  return false;
544 }
545 
546 
548 (
549  const bool report,
550  const scalar minPyrVol,
551  const polyMesh& mesh,
552  const vectorField& cellCentres,
553  const pointField& p,
554  const labelList& checkFaces,
555  const List<labelPair>& baffles,
556  labelHashSet* setPtr
557 )
558 {
559  // check whether face area vector points to the cell with higher label
560  const labelList& own = mesh.faceOwner();
561  const labelList& nei = mesh.faceNeighbour();
562 
563  const faceList& f = mesh.faces();
564 
565  label nErrorPyrs = 0;
566 
567  for (const label facei : checkFaces)
568  {
569  // Create the owner pyramid - it will have negative volume
570  scalar pyrVol = pyramidPointFaceRef
571  (
572  f[facei],
573  cellCentres[own[facei]]
574  ).mag(p);
575 
576  if (pyrVol > -minPyrVol)
577  {
578  ++nErrorPyrs;
579 
580  if (report)
581  {
582  Pout<< "bool polyMeshGeometry::checkFacePyramids("
583  << "const bool, const scalar, const pointField&"
584  << ", const labelList&, labelHashSet*): "
585  << "face " << facei << " points the wrong way. " << endl
586  << "Pyramid volume: " << -pyrVol
587  << " Face " << f[facei] << " area: " << f[facei].mag(p)
588  << " Owner cell: " << own[facei] << endl
589  << "Owner cell vertex labels: "
590  << mesh.cells()[own[facei]].labels(f)
591  << endl;
592  }
593  if (setPtr)
594  {
595  setPtr->insert(facei);
596  }
597  }
598 
599  if (mesh.isInternalFace(facei))
600  {
601  // Create the neighbour pyramid - it will have positive volume
602  scalar pyrVol =
603  pyramidPointFaceRef(f[facei], cellCentres[nei[facei]]).mag(p);
604 
605  if (pyrVol < minPyrVol)
606  {
607  ++nErrorPyrs;
608 
609  if (report)
610  {
611  Pout<< "bool polyMeshGeometry::checkFacePyramids("
612  << "const bool, const scalar, const pointField&"
613  << ", const labelList&, labelHashSet*): "
614  << "face " << facei << " points the wrong way. " << endl
615  << "Pyramid volume: " << -pyrVol
616  << " Face " << f[facei] << " area: " << f[facei].mag(p)
617  << " Neighbour cell: " << nei[facei] << endl
618  << "Neighbour cell vertex labels: "
619  << mesh.cells()[nei[facei]].labels(f)
620  << endl;
621  }
622  if (setPtr)
623  {
624  setPtr->insert(facei);
625  }
626  }
627  }
628  }
629 
630  for (const labelPair& baffle : baffles)
631  {
632  const label face0 = baffle.first();
633  const label face1 = baffle.second();
634 
635  const point& ownCc = cellCentres[own[face0]];
636 
637  // Create the owner pyramid - it will have negative volume
638  scalar pyrVolOwn = pyramidPointFaceRef
639  (
640  f[face0],
641  ownCc
642  ).mag(p);
643 
644  if (pyrVolOwn > -minPyrVol)
645  {
646  ++nErrorPyrs;
647 
648  if (report)
649  {
650  Pout<< "bool polyMeshGeometry::checkFacePyramids("
651  << "const bool, const scalar, const pointField&"
652  << ", const labelList&, labelHashSet*): "
653  << "face " << face0 << " points the wrong way. " << endl
654  << "Pyramid volume: " << -pyrVolOwn
655  << " Face " << f[face0] << " area: " << f[face0].mag(p)
656  << " Owner cell: " << own[face0] << endl
657  << "Owner cell vertex labels: "
658  << mesh.cells()[own[face0]].labels(f)
659  << endl;
660  }
661  if (setPtr)
662  {
663  setPtr->insert(face0);
664  }
665  }
666 
667  // Create the neighbour pyramid - it will have positive volume
668  scalar pyrVolNbr =
669  pyramidPointFaceRef(f[face0], cellCentres[own[face1]]).mag(p);
670 
671  if (pyrVolNbr < minPyrVol)
672  {
673  ++nErrorPyrs;
674 
675  if (report)
676  {
677  Pout<< "bool polyMeshGeometry::checkFacePyramids("
678  << "const bool, const scalar, const pointField&"
679  << ", const labelList&, labelHashSet*): "
680  << "face " << face0 << " points the wrong way. " << endl
681  << "Pyramid volume: " << -pyrVolNbr
682  << " Face " << f[face0] << " area: " << f[face0].mag(p)
683  << " Neighbour cell: " << own[face1] << endl
684  << "Neighbour cell vertex labels: "
685  << mesh.cells()[own[face1]].labels(f)
686  << endl;
687  }
688  if (setPtr)
689  {
690  setPtr->insert(face0);
691  }
692  }
693  }
694 
695  reduce(nErrorPyrs, sumOp<label>());
696 
697  if (nErrorPyrs > 0)
698  {
699  if (report)
700  {
702  << "Error in face pyramids: faces pointing the wrong way."
703  << endl;
704  }
705 
706  return true;
707  }
708 
709  if (report)
710  {
711  Info<< "Face pyramids OK.\n" << endl;
712  }
713 
714  return false;
715 }
716 
717 
719 (
720  const bool report,
721  const scalar minTetQuality,
722  const polyMesh& mesh,
723  const vectorField& cellCentres,
724  const vectorField& faceCentres,
725  const pointField& p,
726  const labelList& checkFaces,
727  const List<labelPair>& baffles,
728  labelHashSet* setPtr
729 )
730 {
731  // check whether decomposing each cell into tets results in
732  // positive volume, non-flat tets
733  const labelList& own = mesh.faceOwner();
734  const labelList& nei = mesh.faceNeighbour();
735  const polyBoundaryMesh& patches = mesh.boundaryMesh();
736 
737  // Calculate coupled cell centre
738  pointField neiCc(mesh.nBoundaryFaces());
739 
740  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
741  {
742  neiCc[facei - mesh.nInternalFaces()] = cellCentres[own[facei]];
743  }
744 
746 
747  label nErrorTets = 0;
748 
749  for (const label facei : checkFaces)
750  {
751  // Create the owner pyramid - note: exchange cell and face centre
752  // to get positive volume.
753  bool tetError = checkFaceTet
754  (
755  mesh,
756  report,
757  minTetQuality,
758  p,
759  facei,
760  cellCentres[own[facei]], // face centre
761  faceCentres[facei], // cell centre
762  setPtr
763  );
764 
765  if (tetError)
766  {
767  ++nErrorTets;
768  }
769 
770  if (mesh.isInternalFace(facei))
771  {
772  // Create the neighbour tets - they will have positive volume
773  bool tetError = checkFaceTet
774  (
775  mesh,
776  report,
777  minTetQuality,
778  p,
779  facei,
780  faceCentres[facei], // face centre
781  cellCentres[nei[facei]], // cell centre
782  setPtr
783  );
784 
785  if (tetError)
786  {
787  ++nErrorTets;
788  }
789 
790  if
791  (
793  (
794  mesh,
795  facei,
796  minTetQuality,
797  report
798  ) == -1
799  )
800  {
801  ++nErrorTets;
802  if (setPtr)
803  {
804  setPtr->insert(facei);
805  }
806  }
807  }
808  else
809  {
810  label patchi = patches.whichPatch(facei);
811 
812  if (patches[patchi].coupled())
813  {
814  if
815  (
817  (
818  mesh,
819  facei,
820  neiCc[facei - mesh.nInternalFaces()],
821  minTetQuality,
822  report
823  ) == -1
824  )
825  {
826  ++nErrorTets;
827  if (setPtr)
828  {
829  setPtr->insert(facei);
830  }
831  }
832  }
833  else
834  {
835  if
836  (
838  (
839  mesh,
840  facei,
841  minTetQuality,
842  report
843  ) == -1
844  )
845  {
846  ++nErrorTets;
847  if (setPtr)
848  {
849  setPtr->insert(facei);
850  }
851  }
852  }
853  }
854  }
855 
856  for (const labelPair& baffle : baffles)
857  {
858  const label face0 = baffle.first();
859  const label face1 = baffle.second();
860 
861  bool tetError = checkFaceTet
862  (
863  mesh,
864  report,
865  minTetQuality,
866  p,
867  face0,
868  cellCentres[own[face0]], // face centre
869  faceCentres[face0], // cell centre
870  setPtr
871  );
872 
873  if (tetError)
874  {
875  ++nErrorTets;
876  }
877 
878  // Create the neighbour tets - they will have positive volume
879  tetError = checkFaceTet
880  (
881  mesh,
882  report,
883  minTetQuality,
884  p,
885  face0,
886  faceCentres[face0], // face centre
887  cellCentres[own[face1]], // cell centre
888  setPtr
889  );
890 
891  if (tetError)
892  {
893  ++nErrorTets;
894  }
895 
896  if
897  (
899  (
900  mesh,
901  face0,
902  cellCentres[own[face1]],
903  minTetQuality,
904  report
905  ) == -1
906  )
907  {
908  ++nErrorTets;
909  if (setPtr)
910  {
911  setPtr->insert(face0);
912  }
913  }
914  }
915 
916  reduce(nErrorTets, sumOp<label>());
917 
918  if (nErrorTets > 0)
919  {
920  if (report)
921  {
923  << "Error in face decomposition: negative tets."
924  << endl;
925  }
926 
927  return true;
928  }
929 
930  if (report)
931  {
932  Info<< "Face tets OK.\n" << endl;
933  }
934 
935  return false;
936 }
937 
938 
940 (
941  const bool report,
942  const scalar internalSkew,
943  const scalar boundarySkew,
944  const polyMesh& mesh,
945  const pointField& points,
946  const vectorField& cellCentres,
947  const vectorField& faceCentres,
948  const vectorField& faceAreas,
949  const labelList& checkFaces,
950  const List<labelPair>& baffles,
951  labelHashSet* setPtr
952 )
953 {
954  // Warn if the skew correction vector is more than skew times
955  // larger than the face area vector
956 
957  const labelList& own = mesh.faceOwner();
958  const labelList& nei = mesh.faceNeighbour();
959  const faceList& faces = mesh.faces();
960  const polyBoundaryMesh& patches = mesh.boundaryMesh();
961 
962  // Calculate coupled cell centre
963  pointField neiCc;
964  syncTools::swapBoundaryCellPositions(mesh, cellCentres, neiCc);
965 
966  scalar maxSkew = 0;
967 
968  label nWarnSkew = 0;
969 
970  for (const label facei : checkFaces)
971  {
972  if (mesh.isInternalFace(facei))
973  {
974  scalar skewness = primitiveMeshTools::faceSkewness
975  (
976  faces,
977  points,
978  faceCentres,
979  faceAreas,
980 
981  facei,
982  cellCentres[own[facei]],
983  cellCentres[nei[facei]]
984  );
985 
986  // Check if the skewness vector is greater than the PN vector.
987  // This does not cause trouble but is a good indication of a poor
988  // mesh.
989  if (skewness > internalSkew)
990  {
991  ++nWarnSkew;
992 
993  if (report)
994  {
995  Pout<< "Severe skewness for face " << facei
996  << " skewness = " << skewness << endl;
997  }
998  if (setPtr)
999  {
1000  setPtr->insert(facei);
1001  }
1002  }
1003 
1004  maxSkew = max(maxSkew, skewness);
1005  }
1006  else if (patches[patches.whichPatch(facei)].coupled())
1007  {
1008  scalar skewness = primitiveMeshTools::faceSkewness
1009  (
1010  faces,
1011  points,
1012  faceCentres,
1013  faceAreas,
1014 
1015  facei,
1016  cellCentres[own[facei]],
1017  neiCc[facei-mesh.nInternalFaces()]
1018  );
1019 
1020  // Check if the skewness vector is greater than the PN vector.
1021  // This does not cause trouble but is a good indication of a poor
1022  // mesh.
1023  if (skewness > internalSkew)
1024  {
1025  ++nWarnSkew;
1026 
1027  if (report)
1028  {
1029  Pout<< "Severe skewness for coupled face " << facei
1030  << " skewness = " << skewness << endl;
1031  }
1032  if (setPtr)
1033  {
1034  setPtr->insert(facei);
1035  }
1036  }
1037 
1038  maxSkew = max(maxSkew, skewness);
1039  }
1040  else
1041  {
1043  (
1044  faces,
1045  points,
1046  faceCentres,
1047  faceAreas,
1048 
1049  facei,
1050  cellCentres[own[facei]]
1051  );
1052 
1053 
1054  // Check if the skewness vector is greater than the PN vector.
1055  // This does not cause trouble but is a good indication of a poor
1056  // mesh.
1057  if (skewness > boundarySkew)
1058  {
1059  ++nWarnSkew;
1060 
1061  if (report)
1062  {
1063  Pout<< "Severe skewness for boundary face " << facei
1064  << " skewness = " << skewness << endl;
1065  }
1066  if (setPtr)
1067  {
1068  setPtr->insert(facei);
1069  }
1070  }
1071 
1072  maxSkew = max(maxSkew, skewness);
1073  }
1074  }
1075 
1076  for (const labelPair& baffle : baffles)
1077  {
1078  const label face0 = baffle.first();
1079  const label face1 = baffle.second();
1080 
1081  const point& ownCc = cellCentres[own[face0]];
1082  const point& neiCc = cellCentres[own[face1]];
1083 
1084  scalar skewness = primitiveMeshTools::faceSkewness
1085  (
1086  faces,
1087  points,
1088  faceCentres,
1089  faceAreas,
1090 
1091  face0,
1092  ownCc,
1093  neiCc
1094  );
1095 
1096  // Check if the skewness vector is greater than the PN vector.
1097  // This does not cause trouble but is a good indication of a poor
1098  // mesh.
1099  if (skewness > internalSkew)
1100  {
1101  ++nWarnSkew;
1102 
1103  if (report)
1104  {
1105  Pout<< "Severe skewness for face " << face0
1106  << " skewness = " << skewness << endl;
1107  }
1108  if (setPtr)
1109  {
1110  setPtr->insert(face0);
1111  }
1112  }
1113 
1114  maxSkew = max(maxSkew, skewness);
1115  }
1116 
1117 
1118  reduce(maxSkew, maxOp<scalar>());
1119  reduce(nWarnSkew, sumOp<label>());
1120 
1121  if (nWarnSkew > 0)
1122  {
1123  if (report)
1124  {
1126  << 100*maxSkew
1127  << " percent.\nThis may impair the quality of the result." << nl
1128  << nWarnSkew << " highly skew faces detected."
1129  << endl;
1130  }
1131 
1132  return true;
1133  }
1134 
1135  if (report)
1136  {
1137  Info<< "Max skewness = " << 100*maxSkew
1138  << " percent. Face skewness OK.\n" << endl;
1139  }
1140 
1141  return false;
1142 }
1143 
1144 
1146 (
1147  const bool report,
1148  const scalar warnWeight,
1149  const polyMesh& mesh,
1150  const vectorField& cellCentres,
1151  const vectorField& faceCentres,
1152  const vectorField& faceAreas,
1153  const labelList& checkFaces,
1154  const List<labelPair>& baffles,
1155  labelHashSet* setPtr
1156 )
1157 {
1158  // Warn if the delta factor (0..1) is too large.
1159 
1160  const labelList& own = mesh.faceOwner();
1161  const labelList& nei = mesh.faceNeighbour();
1162  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1163 
1164  // Calculate coupled cell centre
1165  pointField neiCc(mesh.nBoundaryFaces());
1166 
1167  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
1168  {
1169  neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
1170  }
1172 
1173 
1174  scalar minWeight = GREAT;
1175 
1176  label nWarnWeight = 0;
1177 
1178  for (const label facei : checkFaces)
1179  {
1180  const point& fc = faceCentres[facei];
1181  const vector& fa = faceAreas[facei];
1182 
1183  scalar dOwn = mag(fa & (fc-cellCentres[own[facei]]));
1184 
1185  if (mesh.isInternalFace(facei))
1186  {
1187  scalar dNei = mag(fa & (cellCentres[nei[facei]]-fc));
1188  scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
1189 
1190  if (weight < warnWeight)
1191  {
1192  ++nWarnWeight;
1193 
1194  if (report)
1195  {
1196  Pout<< "Small weighting factor for face " << facei
1197  << " weight = " << weight << endl;
1198  }
1199  if (setPtr)
1200  {
1201  setPtr->insert(facei);
1202  }
1203  }
1204 
1205  minWeight = min(minWeight, weight);
1206  }
1207  else
1208  {
1209  label patchi = patches.whichPatch(facei);
1210 
1211  if (patches[patchi].coupled())
1212  {
1213  scalar dNei = mag(fa & (neiCc[facei-mesh.nInternalFaces()]-fc));
1214  scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
1215 
1216  if (weight < warnWeight)
1217  {
1218  ++nWarnWeight;
1219 
1220  if (report)
1221  {
1222  Pout<< "Small weighting factor for face " << facei
1223  << " weight = " << weight << endl;
1224  }
1225  if (setPtr)
1226  {
1227  setPtr->insert(facei);
1228  }
1229  }
1230 
1231  minWeight = min(minWeight, weight);
1232  }
1233  }
1234  }
1235 
1236  for (const labelPair& baffle : baffles)
1237  {
1238  const label face0 = baffle.first();
1239  const label face1 = baffle.second();
1240 
1241  const point& ownCc = cellCentres[own[face0]];
1242  const point& fc = faceCentres[face0];
1243  const vector& fa = faceAreas[face0];
1244 
1245  scalar dOwn = mag(fa & (fc-ownCc));
1246  scalar dNei = mag(fa & (cellCentres[own[face1]]-fc));
1247  scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
1248 
1249  if (weight < warnWeight)
1250  {
1251  ++nWarnWeight;
1252 
1253  if (report)
1254  {
1255  Pout<< "Small weighting factor for face " << face0
1256  << " weight = " << weight << endl;
1257  }
1258  if (setPtr)
1259  {
1260  setPtr->insert(face0);
1261  }
1262  }
1263 
1264  minWeight = min(minWeight, weight);
1265  }
1266 
1267  reduce(minWeight, minOp<scalar>());
1268  reduce(nWarnWeight, sumOp<label>());
1269 
1270  if (minWeight < warnWeight)
1271  {
1272  if (report)
1273  {
1275  << minWeight << '.' << nl
1276  << nWarnWeight << " faces with small weights detected."
1277  << endl;
1278  }
1279 
1280  return true;
1281  }
1282 
1283  if (report)
1284  {
1285  Info<< "Min weight = " << minWeight
1286  << ". Weights OK.\n" << endl;
1287  }
1288 
1289  return false;
1290 }
1291 
1292 
1294 (
1295  const bool report,
1296  const scalar warnRatio,
1297  const polyMesh& mesh,
1298  const scalarField& cellVolumes,
1299  const labelList& checkFaces,
1300  const List<labelPair>& baffles,
1301  labelHashSet* setPtr
1302 )
1303 {
1304  // Warn if the volume ratio between neighbouring cells is too large
1305 
1306  const labelList& own = mesh.faceOwner();
1307  const labelList& nei = mesh.faceNeighbour();
1308  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1309 
1310  // Calculate coupled cell vol
1311  scalarField neiVols(mesh.nBoundaryFaces());
1312 
1313  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
1314  {
1315  neiVols[facei-mesh.nInternalFaces()] = cellVolumes[own[facei]];
1316  }
1318 
1319 
1320  scalar minRatio = GREAT;
1321 
1322  label nWarnRatio = 0;
1323 
1324  for (const label facei : checkFaces)
1325  {
1326  scalar ownVol = mag(cellVolumes[own[facei]]);
1327 
1328  scalar neiVol = -GREAT;
1329 
1330  if (mesh.isInternalFace(facei))
1331  {
1332  neiVol = mag(cellVolumes[nei[facei]]);
1333  }
1334  else
1335  {
1336  label patchi = patches.whichPatch(facei);
1337 
1338  if (patches[patchi].coupled())
1339  {
1340  neiVol = mag(neiVols[facei-mesh.nInternalFaces()]);
1341  }
1342  }
1343 
1344  if (neiVol >= 0)
1345  {
1346  scalar ratio = min(ownVol, neiVol) / (max(ownVol, neiVol) + VSMALL);
1347 
1348  if (ratio < warnRatio)
1349  {
1350  ++nWarnRatio;
1351 
1352  if (report)
1353  {
1354  Pout<< "Small ratio for face " << facei
1355  << " ratio = " << ratio << endl;
1356  }
1357  if (setPtr)
1358  {
1359  setPtr->insert(facei);
1360  }
1361  }
1362 
1363  minRatio = min(minRatio, ratio);
1364  }
1365  }
1366 
1367  for (const labelPair& baffle : baffles)
1368  {
1369  const label face0 = baffle.first();
1370  const label face1 = baffle.second();
1371 
1372  scalar ownVol = mag(cellVolumes[own[face0]]);
1373 
1374  scalar neiVol = mag(cellVolumes[own[face1]]);
1375 
1376  if (neiVol >= 0)
1377  {
1378  scalar ratio = min(ownVol, neiVol) / (max(ownVol, neiVol) + VSMALL);
1379 
1380  if (ratio < warnRatio)
1381  {
1382  ++nWarnRatio;
1383 
1384  if (report)
1385  {
1386  Pout<< "Small ratio for face " << face0
1387  << " ratio = " << ratio << endl;
1388  }
1389  if (setPtr)
1390  {
1391  setPtr->insert(face0);
1392  }
1393  }
1394 
1395  minRatio = min(minRatio, ratio);
1396  }
1397  }
1398 
1399  reduce(minRatio, minOp<scalar>());
1400  reduce(nWarnRatio, sumOp<label>());
1401 
1402  if (minRatio < warnRatio)
1403  {
1404  if (report)
1405  {
1407  << minRatio << '.' << nl
1408  << nWarnRatio << " faces with small ratios detected."
1409  << endl;
1410  }
1411 
1412  return true;
1413  }
1414 
1415  if (report)
1416  {
1417  Info<< "Min ratio = " << minRatio
1418  << ". Ratios OK.\n" << endl;
1419  }
1420 
1421  return false;
1422 }
1424 
1425 // Check convexity of angles in a face. Allow a slight non-convexity.
1426 // E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity)
1427 // (if truly concave and points not visible from face centre the face-pyramid
1428 // check in checkMesh will fail)
1430 (
1431  const bool report,
1432  const scalar maxDeg,
1433  const polyMesh& mesh,
1434  const vectorField& faceAreas,
1435  const pointField& p,
1436  const labelList& checkFaces,
1437  labelHashSet* setPtr
1438 )
1439 {
1440  if (maxDeg < -SMALL || maxDeg > 180+SMALL)
1441  {
1443  << "maxDeg should be [0..180] but is now " << maxDeg
1444  << abort(FatalError);
1445  }
1446 
1447  const scalar maxSin = Foam::sin(degToRad(maxDeg));
1448 
1449  const faceList& fcs = mesh.faces();
1450 
1451  scalar maxEdgeSin = 0.0;
1452 
1453  label nConcave = 0;
1454 
1455  label errorFacei = -1;
1456 
1457  for (const label facei : checkFaces)
1458  {
1459  const face& f = fcs[facei];
1460 
1461  const vector faceNormal = normalised(faceAreas[facei]);
1462 
1463  // Normalized vector from f[size-1] to f[0];
1464  vector ePrev(p[f.first()] - p[f.last()]);
1465  scalar magEPrev = mag(ePrev);
1466  ePrev /= magEPrev + VSMALL;
1467 
1468  forAll(f, fp0)
1469  {
1470  // Normalized vector between two consecutive points
1471  vector e10(p[f.nextLabel(fp0)] - p[f.thisLabel(fp0)]);
1472  scalar magE10 = mag(e10);
1473  e10 /= magE10 + VSMALL;
1474 
1475  if (magEPrev > SMALL && magE10 > SMALL)
1476  {
1477  vector edgeNormal = ePrev ^ e10;
1478  scalar magEdgeNormal = mag(edgeNormal);
1479 
1480  if (magEdgeNormal < maxSin)
1481  {
1482  // Edges (almost) aligned -> face is ok.
1483  }
1484  else
1485  {
1486  // Check normal
1487  edgeNormal /= magEdgeNormal;
1488 
1489  if ((edgeNormal & faceNormal) < SMALL)
1490  {
1491  if (facei != errorFacei)
1492  {
1493  // Count only one error per face.
1494  errorFacei = facei;
1495  ++nConcave;
1496  }
1497 
1498  maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
1499 
1500  if (setPtr)
1501  {
1502  setPtr->insert(facei);
1503  }
1504  }
1505  }
1506  }
1507 
1508  ePrev = e10;
1509  magEPrev = magE10;
1510  }
1511  }
1512 
1513  reduce(nConcave, sumOp<label>());
1514  reduce(maxEdgeSin, maxOp<scalar>());
1515 
1516  if (report)
1517  {
1518  if (maxEdgeSin > SMALL)
1519  {
1520  scalar maxConcaveDegr =
1521  radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
1522 
1523  Info<< "There are " << nConcave
1524  << " faces with concave angles between consecutive"
1525  << " edges. Max concave angle = " << maxConcaveDegr
1526  << " degrees.\n" << endl;
1527  }
1528  else
1529  {
1530  Info<< "All angles in faces are convex or less than " << maxDeg
1531  << " degrees concave.\n" << endl;
1532  }
1533  }
1534 
1535  if (nConcave > 0)
1536  {
1537  if (report)
1538  {
1540  << nConcave << " face points with severe concave angle (> "
1541  << maxDeg << " deg) found.\n"
1542  << endl;
1543  }
1544 
1545  return true;
1546  }
1547 
1548  return false;
1549 }
1550 
1551 
1552 // Check twist of faces. Is calculated as the difference between normals of
1553 // individual triangles and the cell-cell centre edge.
1555 (
1556  const bool report,
1557  const scalar minTwist,
1558  const polyMesh& mesh,
1559  const vectorField& cellCentres,
1560  const vectorField& faceAreas,
1561  const vectorField& faceCentres,
1562  const pointField& p,
1563  const labelList& checkFaces,
1564  labelHashSet* setPtr
1565 )
1566 {
1567  if (minTwist < -1-SMALL || minTwist > 1+SMALL)
1568  {
1570  << "minTwist should be [-1..1] but is now " << minTwist
1571  << abort(FatalError);
1572  }
1573 
1574 
1575  const faceList& fcs = mesh.faces();
1576 
1577  label nWarped = 0;
1578 
1579  const labelList& own = mesh.faceOwner();
1580  const labelList& nei = mesh.faceNeighbour();
1581  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1582 
1583  // Calculate coupled cell centre
1584  pointField neiCc(mesh.nBoundaryFaces());
1585 
1586  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
1587  {
1588  neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
1589  }
1591 
1592  for (const label facei : checkFaces)
1593  {
1594  const face& f = fcs[facei];
1595 
1596  if (f.size() > 3)
1597  {
1598  vector nf(Zero);
1599 
1600  if (mesh.isInternalFace(facei))
1601  {
1602  nf =
1603  normalised
1604  (
1605  cellCentres[nei[facei]]
1606  - cellCentres[own[facei]]
1607  );
1608  }
1609  else if (patches[patches.whichPatch(facei)].coupled())
1610  {
1611  nf =
1612  normalised
1613  (
1614  neiCc[facei-mesh.nInternalFaces()]
1615  - cellCentres[own[facei]]
1616  );
1617  }
1618  else
1619  {
1620  nf =
1621  normalised
1622  (
1623  faceCentres[facei]
1624  - cellCentres[own[facei]]
1625  );
1626  }
1627 
1628  if (nf != vector::zero)
1629  {
1630  const point& fc = faceCentres[facei];
1631 
1632  forAll(f, fpI)
1633  {
1634  vector triArea
1635  (
1637  (
1638  p[f[fpI]],
1639  p[f.nextLabel(fpI)],
1640  fc
1641  )
1642  );
1643 
1644  scalar magTri = mag(triArea);
1645 
1646  if (magTri > VSMALL && ((nf & triArea/magTri) < minTwist))
1647  {
1648  ++nWarped;
1649 
1650  if (setPtr)
1651  {
1652  setPtr->insert(facei);
1653  }
1654 
1655  break;
1656  }
1657  }
1658  }
1659  }
1660  }
1661 
1662  reduce(nWarped, sumOp<label>());
1663 
1664  if (report)
1665  {
1666  if (nWarped> 0)
1667  {
1668  Info<< "There are " << nWarped
1669  << " faces with cosine of the angle"
1670  << " between triangle normal and face normal less than "
1671  << minTwist << nl << endl;
1672  }
1673  else
1674  {
1675  Info<< "All faces are flat in that the cosine of the angle"
1676  << " between triangle normal and face normal less than "
1677  << minTwist << nl << endl;
1678  }
1679  }
1680 
1681  if (nWarped > 0)
1682  {
1683  if (report)
1684  {
1686  << nWarped << " faces with severe warpage "
1687  << "(cosine of the angle between triangle normal and "
1688  << "face normal < " << minTwist << ") found.\n"
1689  << endl;
1690  }
1691 
1692  return true;
1693  }
1695  return false;
1696 }
1697 
1698 
1699 // Like checkFaceTwist but compares normals of consecutive triangles.
1701 (
1702  const bool report,
1703  const scalar minTwist,
1704  const polyMesh& mesh,
1705  const vectorField& faceAreas,
1706  const vectorField& faceCentres,
1707  const pointField& p,
1708  const labelList& checkFaces,
1709  labelHashSet* setPtr
1710 )
1711 {
1712  if (minTwist < -1-SMALL || minTwist > 1+SMALL)
1713  {
1715  << "minTwist should be [-1..1] but is now " << minTwist
1716  << abort(FatalError);
1717  }
1718 
1719  const faceList& fcs = mesh.faces();
1720 
1721  label nWarped = 0;
1722 
1723  for (const label facei : checkFaces)
1724  {
1725  const face& f = fcs[facei];
1726 
1727  if (f.size() > 3)
1728  {
1729  const point& fc = faceCentres[facei];
1730 
1731  // Find starting triangle (at startFp) with non-zero area
1732  label startFp = -1;
1733  vector prevN;
1734 
1735  forAll(f, fp)
1736  {
1737  prevN = triPointRef::areaNormal
1738  (
1739  p[f[fp]],
1740  p[f.nextLabel(fp)],
1741  fc
1742  );
1743 
1744  scalar magTri = mag(prevN);
1745 
1746  if (magTri > VSMALL)
1747  {
1748  startFp = fp;
1749  prevN /= magTri;
1750  break;
1751  }
1752  }
1753 
1754  if (startFp != -1)
1755  {
1756  label fp = startFp;
1757 
1758  do
1759  {
1760  fp = f.fcIndex(fp);
1761 
1762  vector triN
1763  (
1765  (
1766  p[f[fp]],
1767  p[f.nextLabel(fp)],
1768  fc
1769  )
1770  );
1771  scalar magTri = mag(triN);
1772 
1773  if (magTri > VSMALL)
1774  {
1775  triN /= magTri;
1776 
1777  if ((prevN & triN) < minTwist)
1778  {
1779  ++nWarped;
1780 
1781  if (setPtr)
1782  {
1783  setPtr->insert(facei);
1784  }
1785 
1786  break;
1787  }
1788 
1789  prevN = triN;
1790  }
1791  else if (minTwist > 0)
1792  {
1793  ++nWarped;
1794 
1795  if (setPtr)
1796  {
1797  setPtr->insert(facei);
1798  }
1799 
1800  break;
1801  }
1802  }
1803  while (fp != startFp);
1804  }
1805  }
1806  }
1807 
1808 
1809  reduce(nWarped, sumOp<label>());
1810 
1811  if (report)
1812  {
1813  if (nWarped> 0)
1814  {
1815  Info<< "There are " << nWarped
1816  << " faces with cosine of the angle"
1817  << " between consecutive triangle normals less than "
1818  << minTwist << nl << endl;
1819  }
1820  else
1821  {
1822  Info<< "All faces are flat in that the cosine of the angle"
1823  << " between consecutive triangle normals is less than "
1824  << minTwist << nl << endl;
1825  }
1826  }
1827 
1828  if (nWarped > 0)
1829  {
1830  if (report)
1831  {
1833  << nWarped << " faces with severe warpage "
1834  << "(cosine of the angle between consecutive triangle normals"
1835  << " < " << minTwist << ") found.\n"
1836  << endl;
1837  }
1838 
1839  return true;
1840  }
1841 
1842  return false;
1843 }
1844 
1845 
1847 (
1848  const bool report,
1849  const scalar minFlatness,
1850  const polyMesh& mesh,
1851  const vectorField& faceAreas,
1852  const vectorField& faceCentres,
1853  const pointField& p,
1854  const labelList& checkFaces,
1855  labelHashSet* setPtr
1856 )
1857 {
1858  if (minFlatness < -SMALL || minFlatness > 1+SMALL)
1859  {
1861  << "minFlatness should be [0..1] but is now " << minFlatness
1862  << abort(FatalError);
1863  }
1864 
1865  const faceList& fcs = mesh.faces();
1866 
1867  label nWarped = 0;
1868 
1869  for (const label facei : checkFaces)
1870  {
1871  const face& f = fcs[facei];
1872 
1873  if (f.size() > 3)
1874  {
1875  const point& fc = faceCentres[facei];
1876 
1877  // Sum triangle areas
1878  scalar sumArea = 0.0;
1879 
1880  forAll(f, fp)
1881  {
1882  sumArea += triPointRef
1883  (
1884  p[f[fp]],
1885  p[f.nextLabel(fp)],
1886  fc
1887  ).mag();
1888  }
1889 
1890  if (sumArea/mag(faceAreas[facei]) < minFlatness)
1891  {
1892  ++nWarped;
1893  if (setPtr)
1894  {
1895  setPtr->insert(facei);
1896  }
1897  }
1898  }
1899  }
1900 
1901  reduce(nWarped, sumOp<label>());
1902 
1903  if (report)
1904  {
1905  if (nWarped> 0)
1906  {
1907  Info<< "There are " << nWarped
1908  << " faces with area of individual triangles"
1909  << " compared to overall area less than "
1910  << minFlatness << nl << endl;
1911  }
1912  else
1913  {
1914  Info<< "All faces are flat in that the area of individual triangles"
1915  << " compared to overall area is less than "
1916  << minFlatness << nl << endl;
1917  }
1918  }
1919 
1920  if (nWarped > 0)
1921  {
1922  if (report)
1923  {
1925  << nWarped << " non-flat faces "
1926  << "(area of individual triangles"
1927  << " compared to overall area"
1928  << " < " << minFlatness << ") found.\n"
1929  << endl;
1930  }
1931 
1932  return true;
1933  }
1934 
1935  return false;
1936 }
1937 
1938 
1940 (
1941  const bool report,
1942  const scalar minArea,
1943  const polyMesh& mesh,
1944  const vectorField& faceAreas,
1945  const labelList& checkFaces,
1946  labelHashSet* setPtr
1947 )
1948 {
1949  label nZeroArea = 0;
1950 
1951  for (const label facei : checkFaces)
1952  {
1953  if (mag(faceAreas[facei]) < minArea)
1954  {
1955  ++nZeroArea;
1956  if (setPtr)
1957  {
1958  setPtr->insert(facei);
1959  }
1960  }
1961  }
1962 
1963 
1964  reduce(nZeroArea, sumOp<label>());
1965 
1966  if (report)
1967  {
1968  if (nZeroArea > 0)
1969  {
1970  Info<< "There are " << nZeroArea
1971  << " faces with area < " << minArea << '.' << nl << endl;
1972  }
1973  else
1974  {
1975  Info<< "All faces have area > " << minArea << '.' << nl << endl;
1976  }
1977  }
1978 
1979  if (nZeroArea > 0)
1980  {
1981  if (report)
1982  {
1984  << nZeroArea << " faces with area < " << minArea
1985  << " found.\n"
1986  << endl;
1987  }
1988 
1989  return true;
1990  }
1991 
1992  return false;
1993 }
1994 
1995 
1997 (
1998  const bool report,
1999  const scalar warnDet,
2000  const polyMesh& mesh,
2001  const vectorField& faceAreas,
2002  const labelList& checkFaces,
2003  const labelList& affectedCells,
2004  labelHashSet* setPtr
2005 )
2006 {
2007  const cellList& cells = mesh.cells();
2008 
2009  scalar minDet = GREAT;
2010  scalar sumDet = 0.0;
2011  label nSumDet = 0;
2012  label nWarnDet = 0;
2013 
2014  for (const label celli : affectedCells)
2015  {
2016  const cell& cFaces = cells[celli];
2017 
2018  tensor areaSum(Zero);
2019  scalar magAreaSum = 0;
2020 
2021  for (const label facei : cFaces)
2022  {
2023  scalar magArea = mag(faceAreas[facei]);
2024 
2025  magAreaSum += magArea;
2026  areaSum += faceAreas[facei]*(faceAreas[facei]/(magArea+VSMALL));
2027  }
2028 
2029  scalar scaledDet = det(areaSum/(magAreaSum+VSMALL))/0.037037037037037;
2030 
2031  minDet = min(minDet, scaledDet);
2032  sumDet += scaledDet;
2033  ++nSumDet;
2034 
2035  if (scaledDet < warnDet)
2036  {
2037  ++nWarnDet;
2038  if (setPtr)
2039  {
2040  setPtr->insert(cFaces); // Insert all faces of the cell
2041  }
2042  }
2043  }
2044 
2045  reduce(minDet, minOp<scalar>());
2046  reduce(sumDet, sumOp<scalar>());
2047  reduce(nSumDet, sumOp<label>());
2048  reduce(nWarnDet, sumOp<label>());
2049 
2050  if (report)
2051  {
2052  if (nSumDet > 0)
2053  {
2054  Info<< "Cell determinant (1 = uniform cube) : average = "
2055  << sumDet / nSumDet << " min = " << minDet << endl;
2056  }
2057 
2058  if (nWarnDet > 0)
2059  {
2060  Info<< "There are " << nWarnDet
2061  << " cells with determinant < " << warnDet << '.' << nl
2062  << endl;
2063  }
2064  else
2065  {
2066  Info<< "All faces have determinant > " << warnDet << '.' << nl
2067  << endl;
2068  }
2069  }
2070 
2071  if (nWarnDet > 0)
2072  {
2073  if (report)
2074  {
2076  << nWarnDet << " cells with determinant < " << warnDet
2077  << " found.\n"
2078  << endl;
2079  }
2080 
2081  return true;
2082  }
2083 
2084  return false;
2085 }
2086 
2087 
2089 (
2090  const bool report,
2091  const scalar minEdgeLength,
2092  const polyMesh& mesh,
2093  const labelList& checkFaces,
2094  labelHashSet* setPtr
2095 )
2096 {
2097  const scalar reportLenSqr(Foam::sqr(minEdgeLength));
2098 
2099  const pointField& points = mesh.points();
2100  const faceList& faces = mesh.faces();
2101 
2102  scalar minLenSqr = sqr(GREAT);
2103  scalar maxLenSqr = -sqr(GREAT);
2104 
2105  label nSmall = 0;
2106 
2107  for (const label facei : checkFaces)
2108  {
2109  const face& f = faces[facei];
2110 
2111  forAll(f, fp)
2112  {
2113  label fp1 = f.fcIndex(fp);
2114 
2115  scalar magSqrE = magSqr(points[f[fp]] - points[f[fp1]]);
2116 
2117  if (setPtr && magSqrE < reportLenSqr)
2118  {
2119  if (setPtr->insert(f[fp]) || setPtr->insert(f[fp1]))
2120  {
2121  nSmall++;
2122  }
2123  }
2124 
2125  minLenSqr = min(minLenSqr, magSqrE);
2126  maxLenSqr = max(maxLenSqr, magSqrE);
2127  }
2128  }
2129 
2130  reduce(minLenSqr, minOp<scalar>());
2131  reduce(maxLenSqr, maxOp<scalar>());
2132 
2133  reduce(nSmall, sumOp<label>());
2134 
2135  if (nSmall > 0)
2136  {
2137  if (report)
2138  {
2139  Info<< " *Edges too small, min/max edge length = "
2140  << sqrt(minLenSqr) << " " << sqrt(maxLenSqr)
2141  << ", number too small: " << nSmall << endl;
2142  }
2143 
2144  return true;
2145  }
2146 
2147  if (report)
2148  {
2149  Info<< " Min/max edge length = "
2150  << sqrt(minLenSqr) << " " << sqrt(maxLenSqr) << " OK." << endl;
2151  }
2152 
2153  return false;
2154 }
2155 
2156 
2158 (
2159  const bool report,
2160  const scalar orthWarn,
2161  const labelList& checkFaces,
2162  const List<labelPair>& baffles,
2163  labelHashSet* setPtr
2164 ) const
2165 {
2166  return checkFaceDotProduct
2167  (
2168  report,
2169  orthWarn,
2170  mesh_,
2171  cellCentres_,
2172  faceAreas_,
2173  checkFaces,
2174  baffles,
2175  setPtr
2176  );
2177 }
2178 
2179 
2181 (
2182  const bool report,
2183  const scalar minPyrVol,
2184  const pointField& p,
2185  const labelList& checkFaces,
2186  const List<labelPair>& baffles,
2187  labelHashSet* setPtr
2188 ) const
2189 {
2190  return checkFacePyramids
2191  (
2192  report,
2193  minPyrVol,
2194  mesh_,
2195  cellCentres_,
2196  p,
2197  checkFaces,
2198  baffles,
2199  setPtr
2200  );
2201 }
2202 
2203 
2205 (
2206  const bool report,
2207  const scalar minTetQuality,
2208  const pointField& p,
2209  const labelList& checkFaces,
2210  const List<labelPair>& baffles,
2211  labelHashSet* setPtr
2212 ) const
2213 {
2214  return checkFaceTets
2215  (
2216  report,
2217  minTetQuality,
2218  mesh_,
2219  cellCentres_,
2220  faceCentres_,
2221  p,
2222  checkFaces,
2223  baffles,
2224  setPtr
2225  );
2226 }
2227 
2228 
2229 //bool Foam::polyMeshGeometry::checkFaceSkewness
2230 //(
2231 // const bool report,
2232 // const scalar internalSkew,
2233 // const scalar boundarySkew,
2234 // const labelList& checkFaces,
2235 // const List<labelPair>& baffles,
2236 // labelHashSet* setPtr
2237 //) const
2238 //{
2239 // return checkFaceSkewness
2240 // (
2241 // report,
2242 // internalSkew,
2243 // boundarySkew,
2244 // mesh_,
2245 // cellCentres_,
2246 // faceCentres_,
2247 // faceAreas_,
2248 // checkFaces,
2249 // baffles,
2250 // setPtr
2251 // );
2252 //}
2253 
2254 
2256 (
2257  const bool report,
2258  const scalar warnWeight,
2259  const labelList& checkFaces,
2260  const List<labelPair>& baffles,
2261  labelHashSet* setPtr
2262 ) const
2263 {
2264  return checkFaceWeights
2265  (
2266  report,
2267  warnWeight,
2268  mesh_,
2269  cellCentres_,
2270  faceCentres_,
2271  faceAreas_,
2272  checkFaces,
2273  baffles,
2274  setPtr
2275  );
2276 }
2277 
2278 
2280 (
2281  const bool report,
2282  const scalar warnRatio,
2283  const labelList& checkFaces,
2284  const List<labelPair>& baffles,
2285  labelHashSet* setPtr
2286 ) const
2287 {
2288  return checkVolRatio
2289  (
2290  report,
2291  warnRatio,
2292  mesh_,
2293  cellVolumes_,
2294  checkFaces,
2295  baffles,
2296  setPtr
2297  );
2298 }
2299 
2300 
2302 (
2303  const bool report,
2304  const scalar maxDeg,
2305  const pointField& p,
2306  const labelList& checkFaces,
2307  labelHashSet* setPtr
2308 ) const
2309 {
2310  return checkFaceAngles
2311  (
2312  report,
2313  maxDeg,
2314  mesh_,
2315  faceAreas_,
2316  p,
2317  checkFaces,
2318  setPtr
2319  );
2320 }
2321 
2322 
2324 (
2325  const bool report,
2326  const scalar minTwist,
2327  const pointField& p,
2328  const labelList& checkFaces,
2329  labelHashSet* setPtr
2330 ) const
2331 {
2332  return checkFaceTwist
2333  (
2334  report,
2335  minTwist,
2336  mesh_,
2337  cellCentres_,
2338  faceAreas_,
2339  faceCentres_,
2340  p,
2341  checkFaces,
2342  setPtr
2343  );
2344 }
2345 
2346 
2348 (
2349  const bool report,
2350  const scalar minTwist,
2351  const pointField& p,
2352  const labelList& checkFaces,
2353  labelHashSet* setPtr
2354 ) const
2355 {
2356  return checkTriangleTwist
2357  (
2358  report,
2359  minTwist,
2360  mesh_,
2361  faceAreas_,
2362  faceCentres_,
2363  p,
2364  checkFaces,
2365  setPtr
2366  );
2367 }
2368 
2369 
2371 (
2372  const bool report,
2373  const scalar minFlatness,
2374  const pointField& p,
2375  const labelList& checkFaces,
2376  labelHashSet* setPtr
2377 ) const
2378 {
2379  return checkFaceFlatness
2380  (
2381  report,
2382  minFlatness,
2383  mesh_,
2384  faceAreas_,
2385  faceCentres_,
2386  p,
2387  checkFaces,
2388  setPtr
2389  );
2390 }
2391 
2392 
2394 (
2395  const bool report,
2396  const scalar minArea,
2397  const labelList& checkFaces,
2398  labelHashSet* setPtr
2399 ) const
2400 {
2401  return checkFaceArea
2402  (
2403  report,
2404  minArea,
2405  mesh_,
2406  faceAreas_,
2407  checkFaces,
2408  setPtr
2409  );
2410 }
2411 
2412 
2414 (
2415  const bool report,
2416  const scalar warnDet,
2417  const labelList& checkFaces,
2418  const labelList& affectedCells,
2419  labelHashSet* setPtr
2420 ) const
2421 {
2422  return checkCellDeterminant
2423  (
2424  report,
2425  warnDet,
2426  mesh_,
2427  faceAreas_,
2428  checkFaces,
2429  affectedCells,
2430  setPtr
2431  );
2432 }
2433 
2434 
2436 (
2437  const bool report,
2438  const scalar minEdgeLength,
2439  const labelList& checkFaces,
2440  labelHashSet* setPtr
2441 ) const
2442 {
2443  return checkEdgeLength
2444  (
2445  report,
2446  minEdgeLength,
2447  mesh_,
2448  checkFaces,
2449  setPtr
2450  );
2451 }
2452 
2453 
2454 // ************************************************************************* //
static bool checkFacePyramids(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
See primitiveMesh.
static bool checkFaceTets(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
See primitiveMesh.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensionedScalar acos(const dimensionedScalar &ds)
static bool checkTriangleTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Consecutive triangle (from face-centre decomposition) normals.
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
static bool checkFaceTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Triangle (from face-centre decomposition) normal v.s.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
dimensionedScalar sqrt(const dimensionedScalar &ds)
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const cellList & cells() const
dimensionedScalar det(const dimensionedSphericalTensor &dt)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
static void swapBoundaryCellPositions(const polyMesh &mesh, const UList< point > &cellData, List< point > &neighbourCellData)
Swap to obtain neighbour cell positions for all boundary faces.
Definition: syncTools.C:27
label nFaces() const noexcept
Number of mesh faces.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
static bool checkVolRatio(const bool report, const scalar warnRatio, const polyMesh &mesh, const scalarField &cellVolumes, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
Cell volume ratio of neighbouring cells (1 for regular mesh)
dimensionedScalar asin(const dimensionedScalar &ds)
static bool checkFaceFlatness(const bool report, const scalar minFlatness, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Area of faces v.s. sum of triangle areas.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
static label findBasePoint(const polyMesh &mesh, label fI, scalar tol, bool report=false)
Find the base point to use for a minimum triangle.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:97
pyramid< point, const point &, const face & > pyramidPointFaceRef
A pyramid using referred point and face.
Definition: pyramid.H:70
Vector< solveScalar > solveVector
Definition: vector.H:60
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
static tmp< scalarField > faceSkewness(const primitiveMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
const cellShapeList & cells
const pointField & points
static bool checkEdgeLength(const bool report, const scalar minEdgeLength, const polyMesh &mesh, const labelList &checkFaces, labelHashSet *pointSetPtr)
Small edges. Optionally collects points of small edges.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
Point centre() const
Return centre (centroid)
Definition: triangleI.H:126
constexpr scalar pi(M_PI)
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static bool checkFaceArea(const bool report, const scalar minArea, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Small faces.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:867
static bool checkFaceAngles(const bool report, const scalar maxDeg, const polyMesh &mesh, const vectorField &faceAreas, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
See primitiveMesh.
static scalar boundaryFaceSkewness(const UList< face > &faces, const pointField &p, const vectorField &fCtrs, const vectorField &fAreas, const label facei, const point &ownCc)
Skewness of single boundary face.
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
static bool checkFaceDotProduct(const bool report, const scalar orthWarn, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
See primitiveMesh.
const dimensionedScalar c
Speed of light in a vacuum.
const polyBoundaryMesh & patches
static void swapBoundaryFacePositions(const polyMesh &mesh, UList< point > &positions)
Swap coupled positions. Uses eqOp.
Definition: syncTools.H:503
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
void correct()
Take over properties from mesh.
static labelList affectedCells(const polyMesh &, const labelList &changedFaces)
Helper function: get affected cells from faces.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
static bool checkFaceWeights(const bool report, const scalar warnWeight, const polyMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
Interpolation weights (0.5 for regular mesh)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static bool checkFaceSkewness(const bool report, const scalar internalSkew, const scalar boundarySkew, const polyMesh &mesh, const pointField &points, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
See primitiveMesh.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
bool coupled
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
static bool checkCellDeterminant(const bool report, const scalar minDet, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, const labelList &affectedCells, labelHashSet *setPtr)
Area of internal faces v.s. boundary faces.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:485
polyMeshGeometry(const polyMesh &)
Construct from mesh.
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:151
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127