primitiveMeshGeometry.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-2016 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 
29 #include "primitiveMeshGeometry.H"
30 #include "pyramidPointFaceRef.H"
31 #include "unitConversion.H"
32 #include "triPointRef.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(primitiveMeshGeometry, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::primitiveMeshGeometry::updateFaceCentresAndAreas
45 (
46  const pointField& p,
47  const labelList& changedFaces
48 )
49 {
50  const faceList& fs = mesh_.faces();
51 
52  forAll(changedFaces, i)
53  {
54  label facei = changedFaces[i];
55 
56  const labelList& f = fs[facei];
57  label nPoints = f.size();
58 
59  // If the face is a triangle, do a direct calculation for efficiency
60  // and to avoid round-off error-related problems
61  if (nPoints == 3)
62  {
63  faceCentres_[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
64  faceAreas_[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
65  }
66  else
67  {
68  vector sumN = Zero;
69  scalar sumA = 0.0;
70  vector sumAc = Zero;
71 
72  point fCentre = p[f[0]];
73  for (label pi = 1; pi < nPoints; pi++)
74  {
75  fCentre += p[f[pi]];
76  }
77 
78  fCentre /= nPoints;
79 
80  for (label pi = 0; pi < nPoints; pi++)
81  {
82  const point& nextPoint = p[f[(pi + 1) % nPoints]];
83 
84  vector c = p[f[pi]] + nextPoint + fCentre;
85  vector n = (nextPoint - p[f[pi]])^(fCentre - p[f[pi]]);
86  scalar a = mag(n);
87 
88  sumN += n;
89  sumA += a;
90  sumAc += a*c;
91  }
92 
93  faceCentres_[facei] = (1.0/3.0)*sumAc/(sumA + VSMALL);
94  faceAreas_[facei] = 0.5*sumN;
95  }
96  }
97 }
98 
99 
100 void Foam::primitiveMeshGeometry::updateCellCentresAndVols
101 (
102  const labelList& changedCells,
103  const labelList& changedFaces
104 )
105 {
106  // Clear the fields for accumulation
107  UIndirectList<vector>(cellCentres_, changedCells) = Zero;
108  UIndirectList<scalar>(cellVolumes_, changedCells) = 0.0;
109 
110  const labelList& own = mesh_.faceOwner();
111  const labelList& nei = mesh_.faceNeighbour();
112 
113  // first estimate the approximate cell centre as the average of face centres
114 
115  vectorField cEst(mesh_.nCells());
116  UIndirectList<vector>(cEst, changedCells) = Zero;
117  scalarField nCellFaces(mesh_.nCells());
118  UIndirectList<scalar>(nCellFaces, changedCells) = 0.0;
119 
120  forAll(changedFaces, i)
121  {
122  label facei = changedFaces[i];
123  cEst[own[facei]] += faceCentres_[facei];
124  nCellFaces[own[facei]] += 1;
125 
126  if (mesh_.isInternalFace(facei))
127  {
128  cEst[nei[facei]] += faceCentres_[facei];
129  nCellFaces[nei[facei]] += 1;
130  }
131  }
132 
133  forAll(changedCells, i)
134  {
135  label celli = changedCells[i];
136  cEst[celli] /= nCellFaces[celli];
137  }
138 
139  forAll(changedFaces, i)
140  {
141  label facei = changedFaces[i];
142 
143  // Calculate 3*face-pyramid volume
144  scalar pyr3Vol = max
145  (
146  faceAreas_[facei] & (faceCentres_[facei] - cEst[own[facei]]),
147  VSMALL
148  );
149 
150  // Calculate face-pyramid centre
151  vector pc = (3.0/4.0)*faceCentres_[facei] + (1.0/4.0)*cEst[own[facei]];
152 
153  // Accumulate volume-weighted face-pyramid centre
154  cellCentres_[own[facei]] += pyr3Vol*pc;
155 
156  // Accumulate face-pyramid volume
157  cellVolumes_[own[facei]] += pyr3Vol;
158 
159  if (mesh_.isInternalFace(facei))
160  {
161  // Calculate 3*face-pyramid volume
162  scalar pyr3Vol = max
163  (
164  faceAreas_[facei] & (cEst[nei[facei]] - faceCentres_[facei]),
165  VSMALL
166  );
167 
168  // Calculate face-pyramid centre
169  vector pc =
170  (3.0/4.0)*faceCentres_[facei]
171  + (1.0/4.0)*cEst[nei[facei]];
172 
173  // Accumulate volume-weighted face-pyramid centre
174  cellCentres_[nei[facei]] += pyr3Vol*pc;
175 
176  // Accumulate face-pyramid volume
177  cellVolumes_[nei[facei]] += pyr3Vol;
178  }
179  }
180 
181  forAll(changedCells, i)
182  {
183  label celli = changedCells[i];
184 
185  cellCentres_[celli] /= cellVolumes_[celli];
186  cellVolumes_[celli] *= (1.0/3.0);
187  }
188 }
189 
190 
192 (
193  const labelList& changedFaces
194 ) const
195 {
196  const labelList& own = mesh_.faceOwner();
197  const labelList& nei = mesh_.faceNeighbour();
198 
199  labelHashSet affectedCells(2*changedFaces.size());
200 
201  forAll(changedFaces, i)
202  {
203  label facei = changedFaces[i];
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 
217 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
218 
219 // Construct from components
221 (
222  const primitiveMesh& mesh
223 )
224 :
225  mesh_(mesh)
226 {
227  correct();
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 {
238  faceAreas_ = mesh_.faceAreas();
239  faceCentres_ = mesh_.faceCentres();
240  cellCentres_ = mesh_.cellCentres();
241  cellVolumes_ = mesh_.cellVolumes();
242 }
243 
244 
246 (
247  const pointField& p,
248  const labelList& changedFaces
249 )
250 {
251  // Update face quantities
252  updateFaceCentresAndAreas(p, changedFaces);
253  // Update cell quantities from face quantities
254  updateCellCentresAndVols(affectedCells(changedFaces), changedFaces);
255 }
256 
257 
259 (
260  const bool report,
261  const scalar orthWarn,
262  const primitiveMesh& mesh,
263  const vectorField& cellCentres,
264  const vectorField& faceAreas,
265  const labelList& checkFaces,
266  labelHashSet* setPtr
267 )
268 {
269  // for all internal faces check that the d dot S product is positive
270 
271  const labelList& own = mesh.faceOwner();
272  const labelList& nei = mesh.faceNeighbour();
273 
274  // Severe nonorthogonality threshold
275  const scalar severeNonorthogonalityThreshold = ::cos(degToRad(orthWarn));
276 
277  scalar minDDotS = GREAT;
278 
279  scalar sumDDotS = 0;
280 
281  label severeNonOrth = 0;
282 
283  label errorNonOrth = 0;
284 
285  forAll(checkFaces, i)
286  {
287  label facei = checkFaces[i];
288 
289  if (mesh.isInternalFace(facei))
290  {
291  vector d = cellCentres[nei[facei]] - cellCentres[own[facei]];
292  const vector& s = faceAreas[facei];
293 
294  scalar dDotS = (d & s)/(mag(d)*mag(s) + VSMALL);
295 
296  if (dDotS < severeNonorthogonalityThreshold)
297  {
298  if (dDotS > SMALL)
299  {
300  if (report)
301  {
302  // Severe non-orthogonality but mesh still OK
303  Pout<< "Severe non-orthogonality for face " << facei
304  << " between cells " << own[facei]
305  << " and " << nei[facei]
306  << ": Angle = " << radToDeg(::acos(dDotS))
307  << " deg." << endl;
308  }
309 
310  if (setPtr)
311  {
312  setPtr->insert(facei);
313  }
314 
315  severeNonOrth++;
316  }
317  else
318  {
319  // Non-orthogonality greater than 90 deg
320  if (report)
321  {
323  << "Severe non-orthogonality detected for face "
324  << facei
325  << " between cells " << own[facei] << " and "
326  << nei[facei]
327  << ": Angle = " << radToDeg(::acos(dDotS))
328  << " deg." << endl;
329  }
330 
331  errorNonOrth++;
332 
333  if (setPtr)
334  {
335  setPtr->insert(facei);
336  }
337  }
338  }
339 
340  if (dDotS < minDDotS)
341  {
342  minDDotS = dDotS;
343  }
344 
345  sumDDotS += dDotS;
346  }
347  }
348 
349  reduce(minDDotS, minOp<scalar>());
350  reduce(sumDDotS, sumOp<scalar>());
351  reduce(severeNonOrth, sumOp<label>());
352  reduce(errorNonOrth, sumOp<label>());
353 
354  label neiSize = nei.size();
355  reduce(neiSize, sumOp<label>());
356 
357  // Only report if there are some internal faces
358  if (neiSize > 0)
359  {
360  if (report && minDDotS < severeNonorthogonalityThreshold)
361  {
362  Info<< "Number of non-orthogonality errors: " << errorNonOrth
363  << ". Number of severely non-orthogonal faces: "
364  << severeNonOrth << "." << endl;
365  }
366  }
367 
368  if (report)
369  {
370  if (neiSize > 0)
371  {
372  Info<< "Mesh non-orthogonality Max: "
373  << radToDeg(::acos(minDDotS))
374  << " average: " << radToDeg(::acos(sumDDotS/neiSize))
375  << endl;
376  }
377  }
378 
379  if (errorNonOrth > 0)
380  {
381  if (report)
382  {
384  << "Error in non-orthogonality detected" << endl;
385  }
386 
387  return true;
388  }
389 
390  if (report)
391  {
392  Info<< "Non-orthogonality check OK.\n" << endl;
393  }
394 
395  return false;
396 }
397 
398 
400 (
401  const bool report,
402  const scalar minPyrVol,
403  const primitiveMesh& mesh,
404  const vectorField& cellCentres,
405  const pointField& p,
406  const labelList& checkFaces,
407  labelHashSet* setPtr
408 )
409 {
410  // check whether face area vector points to the cell with higher label
411  const labelList& own = mesh.faceOwner();
412  const labelList& nei = mesh.faceNeighbour();
413 
414  const faceList& f = mesh.faces();
415 
416  label nErrorPyrs = 0;
417 
418  forAll(checkFaces, i)
419  {
420  label facei = checkFaces[i];
421 
422  // Create the owner pyramid - it will have negative volume
423  scalar pyrVol = pyramidPointFaceRef
424  (
425  f[facei],
426  cellCentres[own[facei]]
427  ).mag(p);
428 
429  if (pyrVol > -minPyrVol)
430  {
431  if (report)
432  {
433  Pout<< "bool primitiveMeshGeometry::checkFacePyramids("
434  << "const bool, const scalar, const pointField&"
435  << ", const labelList&, labelHashSet*): "
436  << "face " << facei << " points the wrong way. " << endl
437  << "Pyramid volume: " << -pyrVol
438  << " Face " << f[facei] << " area: " << f[facei].mag(p)
439  << " Owner cell: " << own[facei] << endl
440  << "Owner cell vertex labels: "
441  << mesh.cells()[own[facei]].labels(f)
442  << endl;
443  }
444 
445 
446  if (setPtr)
447  {
448  setPtr->insert(facei);
449  }
450 
451  nErrorPyrs++;
452  }
453 
454  if (mesh.isInternalFace(facei))
455  {
456  // Create the neighbour pyramid - it will have positive volume
457  scalar pyrVol =
458  pyramidPointFaceRef(f[facei], cellCentres[nei[facei]]).mag(p);
459 
460  if (pyrVol < minPyrVol)
461  {
462  if (report)
463  {
464  Pout<< "bool primitiveMeshGeometry::checkFacePyramids("
465  << "const bool, const scalar, const pointField&"
466  << ", const labelList&, labelHashSet*): "
467  << "face " << facei << " points the wrong way. " << endl
468  << "Pyramid volume: " << -pyrVol
469  << " Face " << f[facei] << " area: " << f[facei].mag(p)
470  << " Neighbour cell: " << nei[facei] << endl
471  << "Neighbour cell vertex labels: "
472  << mesh.cells()[nei[facei]].labels(f)
473  << endl;
474  }
475 
476  if (setPtr)
477  {
478  setPtr->insert(facei);
479  }
480 
481  nErrorPyrs++;
482  }
483  }
484  }
485 
486  reduce(nErrorPyrs, sumOp<label>());
487 
488  if (nErrorPyrs > 0)
489  {
490  if (report)
491  {
493  << "Error in face pyramids: faces pointing the wrong way!"
494  << endl;
495  }
496 
497  return true;
498  }
499 
500  if (report)
501  {
502  Info<< "Face pyramids OK.\n" << endl;
503  }
504 
505  return false;
506 }
507 
508 
510 (
511  const bool report,
512  const scalar internalSkew,
513  const scalar boundarySkew,
514  const primitiveMesh& mesh,
515  const vectorField& cellCentres,
516  const vectorField& faceCentres,
517  const vectorField& faceAreas,
518  const labelList& checkFaces,
519  labelHashSet* setPtr
520 )
521 {
522  // Warn if the skew correction vector is more than skew times
523  // larger than the face area vector
524 
525  const labelList& own = mesh.faceOwner();
526  const labelList& nei = mesh.faceNeighbour();
527 
528  scalar maxSkew = 0;
529 
530  label nWarnSkew = 0;
531 
532  forAll(checkFaces, i)
533  {
534  label facei = checkFaces[i];
535 
536  if (mesh.isInternalFace(facei))
537  {
538  scalar dOwn = mag(faceCentres[facei] - cellCentres[own[facei]]);
539  scalar dNei = mag(faceCentres[facei] - cellCentres[nei[facei]]);
540 
541  point faceIntersection =
542  cellCentres[own[facei]]*dNei/(dOwn+dNei)
543  + cellCentres[nei[facei]]*dOwn/(dOwn+dNei);
544 
545  scalar skewness =
546  mag(faceCentres[facei] - faceIntersection)
547  / (
548  mag(cellCentres[nei[facei]]-cellCentres[own[facei]])
549  + VSMALL
550  );
551 
552  // Check if the skewness vector is greater than the PN vector.
553  // This does not cause trouble but is a good indication of a poor
554  // mesh.
555  if (skewness > internalSkew)
556  {
557  if (report)
558  {
559  Pout<< "Severe skewness for face " << facei
560  << " skewness = " << skewness << endl;
561  }
562 
563  if (setPtr)
564  {
565  setPtr->insert(facei);
566  }
567 
568  nWarnSkew++;
569  }
570 
571  if (skewness > maxSkew)
572  {
573  maxSkew = skewness;
574  }
575  }
576  else
577  {
578  // Boundary faces: consider them to have only skewness error.
579  // (i.e. treat as if mirror cell on other side)
580 
581  const vector faceNormal = normalised(faceAreas[facei]);
582 
583  vector dOwn = faceCentres[facei] - cellCentres[own[facei]];
584 
585  vector dWall = faceNormal*(faceNormal & dOwn);
586 
587  point faceIntersection = cellCentres[own[facei]] + dWall;
588 
589  scalar skewness =
590  mag(faceCentres[facei] - faceIntersection)
591  /(2*mag(dWall) + VSMALL);
592 
593  // Check if the skewness vector is greater than the PN vector.
594  // This does not cause trouble but is a good indication of a poor
595  // mesh.
596  if (skewness > boundarySkew)
597  {
598  if (report)
599  {
600  Pout<< "Severe skewness for boundary face " << facei
601  << " skewness = " << skewness << endl;
602  }
603 
604  if (setPtr)
605  {
606  setPtr->insert(facei);
607  }
608 
609  nWarnSkew++;
610  }
611 
612  if (skewness > maxSkew)
613  {
614  maxSkew = skewness;
615  }
616  }
617  }
618 
619  reduce(maxSkew, maxOp<scalar>());
620  reduce(nWarnSkew, sumOp<label>());
621 
622  if (nWarnSkew > 0)
623  {
624  if (report)
625  {
627  << 100*maxSkew
628  << " percent.\nThis may impair the quality of the result." << nl
629  << nWarnSkew << " highly skew faces detected."
630  << endl;
631  }
632 
633  return true;
634  }
635 
636  if (report)
637  {
638  Info<< "Max skewness = " << 100*maxSkew
639  << " percent. Face skewness OK.\n" << endl;
640  }
641 
642  return false;
643 }
644 
645 
647 (
648  const bool report,
649  const scalar warnWeight,
650  const primitiveMesh& mesh,
651  const vectorField& cellCentres,
652  const vectorField& faceCentres,
653  const vectorField& faceAreas,
654  const labelList& checkFaces,
655  labelHashSet* setPtr
656 )
657 {
658  // Warn if the delta factor (0..1) is too large.
659 
660  const labelList& own = mesh.faceOwner();
661  const labelList& nei = mesh.faceNeighbour();
662 
663  scalar minWeight = GREAT;
664 
665  label nWarnWeight = 0;
666 
667  forAll(checkFaces, i)
668  {
669  label facei = checkFaces[i];
670 
671  if (mesh.isInternalFace(facei))
672  {
673  const point& fc = faceCentres[facei];
674 
675  scalar dOwn = mag(faceAreas[facei] & (fc-cellCentres[own[facei]]));
676  scalar dNei = mag(faceAreas[facei] & (cellCentres[nei[facei]]-fc));
677 
678  scalar weight = min(dNei,dOwn)/(dNei+dOwn);
679 
680  if (weight < warnWeight)
681  {
682  if (report)
683  {
684  Pout<< "Small weighting factor for face " << facei
685  << " weight = " << weight << endl;
686  }
687 
688  if (setPtr)
689  {
690  setPtr->insert(facei);
691  }
692 
693  nWarnWeight++;
694  }
695 
696  minWeight = min(minWeight, weight);
697  }
698  }
699 
700  reduce(minWeight, minOp<scalar>());
701  reduce(nWarnWeight, sumOp<label>());
702 
703  if (minWeight < warnWeight)
704  {
705  if (report)
706  {
708  << minWeight << '.' << nl
709  << nWarnWeight << " faces with small weights detected."
710  << endl;
711  }
712 
713  return true;
714  }
715 
716  if (report)
717  {
718  Info<< "Min weight = " << minWeight
719  << " percent. Weights OK.\n" << endl;
720  }
721 
722  return false;
723 }
724 
725 
726 // Check convexity of angles in a face. Allow a slight non-convexity.
727 // E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity)
728 // (if truly concave and points not visible from face centre the face-pyramid
729 // check in checkMesh will fail)
731 (
732  const bool report,
733  const scalar maxDeg,
734  const primitiveMesh& mesh,
735  const vectorField& faceAreas,
736  const pointField& p,
737  const labelList& checkFaces,
738  labelHashSet* setPtr
739 )
740 {
741  if (maxDeg < -SMALL || maxDeg > 180+SMALL)
742  {
744  << "maxDeg should be [0..180] but is now " << maxDeg
745  << abort(FatalError);
746  }
747 
748  const scalar maxSin = Foam::sin(degToRad(maxDeg));
749 
750  const faceList& fcs = mesh.faces();
751 
752  scalar maxEdgeSin = 0.0;
753 
754  label nConcave = 0;
755 
756  label errorFacei = -1;
757 
758  forAll(checkFaces, i)
759  {
760  label facei = checkFaces[i];
761 
762  const face& f = fcs[facei];
763 
764  const vector faceNormal = normalised(faceAreas[facei]);
765 
766  // Get edge from f[0] to f[size-1];
767  vector ePrev(p[f.first()] - p[f.last()]);
768  scalar magEPrev = mag(ePrev);
769  ePrev /= magEPrev + VSMALL;
770 
771  forAll(f, fp0)
772  {
773  // Get vertex after fp
774  label fp1 = f.fcIndex(fp0);
775 
776  // Normalized vector between two consecutive points
777  vector e10(p[f[fp1]] - p[f[fp0]]);
778  scalar magE10 = mag(e10);
779  e10 /= magE10 + VSMALL;
780 
781  if (magEPrev > SMALL && magE10 > SMALL)
782  {
783  vector edgeNormal = ePrev ^ e10;
784  scalar magEdgeNormal = mag(edgeNormal);
785 
786  if (magEdgeNormal < maxSin)
787  {
788  // Edges (almost) aligned -> face is ok.
789  }
790  else
791  {
792  // Check normal
793  edgeNormal /= magEdgeNormal;
794 
795  if ((edgeNormal & faceNormal) < SMALL)
796  {
797  if (facei != errorFacei)
798  {
799  // Count only one error per face.
800  errorFacei = facei;
801  nConcave++;
802  }
803 
804  if (setPtr)
805  {
806  setPtr->insert(facei);
807  }
808 
809  maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
810  }
811  }
812  }
813 
814  ePrev = e10;
815  magEPrev = magE10;
816  }
817  }
818 
819  reduce(nConcave, sumOp<label>());
820  reduce(maxEdgeSin, maxOp<scalar>());
821 
822  if (report)
823  {
824  if (maxEdgeSin > SMALL)
825  {
826  scalar maxConcaveDegr =
827  radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
828 
829  Info<< "There are " << nConcave
830  << " faces with concave angles between consecutive"
831  << " edges. Max concave angle = " << maxConcaveDegr
832  << " degrees.\n" << endl;
833  }
834  else
835  {
836  Info<< "All angles in faces are convex or less than " << maxDeg
837  << " degrees concave.\n" << endl;
838  }
839  }
840 
841  if (nConcave > 0)
842  {
843  if (report)
844  {
846  << nConcave << " face points with severe concave angle (> "
847  << maxDeg << " deg) found.\n"
848  << endl;
849  }
850 
851  return true;
852  }
853 
854  return false;
855 }
856 
857 
861 //bool Foam::primitiveMeshGeometry::checkFaceFlatness
862 //(
863 // const bool report,
864 // const scalar warnFlatness,
865 // const primitiveMesh& mesh,
866 // const vectorField& faceAreas,
867 // const vectorField& faceCentres,
868 // const pointField& p,
869 // const labelList& checkFaces,
870 // labelHashSet* setPtr
871 //)
872 //{
873 // if (warnFlatness < 0 || warnFlatness > 1)
874 // {
875 // FatalErrorInFunction
876 // << "warnFlatness should be [0..1] but is now " << warnFlatness
877 // << abort(FatalError);
878 // }
879 //
880 //
881 // const faceList& fcs = mesh.faces();
882 //
883 // // Areas are calculated as the sum of areas. (see
884 // // primitiveMeshFaceCentresAndAreas.C)
885 //
886 // label nWarped = 0;
887 //
888 // scalar minFlatness = GREAT;
889 // scalar sumFlatness = 0;
890 // label nSummed = 0;
891 //
892 // forAll(checkFaces, i)
893 // {
894 // label facei = checkFaces[i];
895 //
896 // const face& f = fcs[facei];
897 //
898 // scalar magArea = mag(faceAreas[facei]);
899 //
900 // if (f.size() > 3 && magArea > VSMALL)
901 // {
902 // const point& fc = faceCentres[facei];
903 //
904 // // Calculate the sum of magnitude of areas and compare to
905 // // magnitude of sum of areas.
906 //
907 // scalar sumA = 0.0;
908 //
909 // forAll(f, fp)
910 // {
911 // const point& thisPoint = p[f[fp]];
912 // const point& nextPoint = p[f.nextLabel(fp)];
913 //
914 // // Triangle around fc.
915 // vector n = 0.5*((nextPoint - thisPoint)^(fc - thisPoint));
916 // sumA += mag(n);
917 // }
918 //
919 // scalar flatness = magArea / (sumA+VSMALL);
920 //
921 // sumFlatness += flatness;
922 // nSummed++;
923 //
924 // minFlatness = min(minFlatness, flatness);
925 //
926 // if (flatness < warnFlatness)
927 // {
928 // nWarped++;
929 //
930 // if (setPtr)
931 // {
932 // setPtr->insert(facei);
933 // }
934 // }
935 // }
936 // }
937 //
938 //
939 // reduce(nWarped, sumOp<label>());
940 // reduce(minFlatness, minOp<scalar>());
941 //
942 // reduce(nSummed, sumOp<label>());
943 // reduce(sumFlatness, sumOp<scalar>());
944 //
945 // if (report)
946 // {
947 // if (nSummed > 0)
948 // {
949 // Info<< "Face flatness (1 = flat, 0 = butterfly) : average = "
950 // << sumFlatness / nSummed << " min = " << minFlatness << endl;
951 // }
952 //
953 // if (nWarped> 0)
954 // {
955 // Info<< "There are " << nWarped
956 // << " faces with ratio between projected and actual area < "
957 // << warnFlatness
958 // << ".\nMinimum ratio (minimum flatness, maximum warpage) = "
959 // << minFlatness << nl << endl;
960 // }
961 // else
962 // {
963 // Info<< "All faces are flat in that the ratio between projected"
964 // << " and actual area is > " << warnFlatness << nl << endl;
965 // }
966 // }
967 //
968 // if (nWarped > 0)
969 // {
970 // if (report)
971 // {
972 // WarningInFunction
973 // << nWarped << " faces with severe warpage (flatness < "
974 // << warnFlatness << ") found.\n"
975 // << endl;
976 // }
977 //
978 // return true;
979 // }
980 //
981 // return false;
982 //}
983 
984 
985 // Check twist of faces. Is calculated as the difference between areas of
986 // individual triangles and the overall area of the face (which ifself is
987 // is the average of the areas of the individual triangles).
989 (
990  const bool report,
991  const scalar minTwist,
992  const primitiveMesh& mesh,
993  const vectorField& faceAreas,
994  const vectorField& faceCentres,
995  const pointField& p,
996  const labelList& checkFaces,
997  labelHashSet* setPtr
998 )
999 {
1000  if (minTwist < -1-SMALL || minTwist > 1+SMALL)
1001  {
1003  << "minTwist should be [-1..1] but is now " << minTwist
1004  << abort(FatalError);
1005  }
1006 
1007 
1008  const faceList& fcs = mesh.faces();
1009 
1010  // Areas are calculated as the sum of areas. (see
1011  // primitiveMeshFaceCentresAndAreas.C)
1012 
1013  label nWarped = 0;
1014 
1015  for (const label facei : checkFaces)
1016  {
1017  const face& f = fcs[facei];
1018 
1019  const scalar magArea = mag(faceAreas[facei]);
1020 
1021  if (f.size() > 3 && magArea > VSMALL)
1022  {
1023  const vector nf = faceAreas[facei] / magArea;
1024 
1025  const point& fc = faceCentres[facei];
1026 
1027  forAll(f, fpI)
1028  {
1029  const vector triArea
1030  (
1031  triPointRef
1032  (
1033  p[f[fpI]],
1034  p[f.nextLabel(fpI)],
1035  fc
1036  ).areaNormal()
1037  );
1038 
1039  const scalar magTri = mag(triArea);
1040 
1041  if (magTri > VSMALL && ((nf & triArea/magTri) < minTwist))
1042  {
1043  ++nWarped;
1044 
1045  if (setPtr)
1046  {
1047  setPtr->insert(facei);
1048  }
1049  }
1050  }
1051  }
1052  }
1053 
1054 
1055  reduce(nWarped, sumOp<label>());
1056 
1057  if (report)
1058  {
1059  if (nWarped> 0)
1060  {
1061  Info<< "There are " << nWarped
1062  << " faces with cosine of the angle"
1063  << " between triangle normal and face normal less than "
1064  << minTwist << nl << endl;
1065  }
1066  else
1067  {
1068  Info<< "All faces are flat in that the cosine of the angle"
1069  << " between triangle normal and face normal less than "
1070  << minTwist << nl << endl;
1071  }
1072  }
1073 
1074  if (nWarped > 0)
1075  {
1076  if (report)
1077  {
1079  << nWarped << " faces with severe warpage "
1080  << "(cosine of the angle between triangle normal and "
1081  << "face normal < " << minTwist << ") found.\n"
1082  << endl;
1083  }
1084 
1085  return true;
1086  }
1087 
1088  return false;
1089 }
1090 
1091 
1094  const bool report,
1095  const scalar minArea,
1096  const primitiveMesh& mesh,
1097  const vectorField& faceAreas,
1098  const labelList& checkFaces,
1099  labelHashSet* setPtr
1100 )
1101 {
1102  label nZeroArea = 0;
1103 
1104  forAll(checkFaces, i)
1105  {
1106  label facei = checkFaces[i];
1107 
1108  if (mag(faceAreas[facei]) < minArea)
1109  {
1110  if (setPtr)
1111  {
1112  setPtr->insert(facei);
1113  }
1114  nZeroArea++;
1115  }
1116  }
1117 
1118 
1119  reduce(nZeroArea, sumOp<label>());
1120 
1121  if (report)
1122  {
1123  if (nZeroArea > 0)
1124  {
1125  Info<< "There are " << nZeroArea
1126  << " faces with area < " << minArea << '.' << nl << endl;
1127  }
1128  else
1129  {
1130  Info<< "All faces have area > " << minArea << '.' << nl << endl;
1131  }
1132  }
1133 
1134  if (nZeroArea > 0)
1135  {
1136  if (report)
1137  {
1139  << nZeroArea << " faces with area < " << minArea
1140  << " found.\n"
1141  << endl;
1142  }
1143 
1144  return true;
1145  }
1146 
1147  return false;
1148 }
1149 
1150 
1153  const bool report,
1154  const scalar warnDet,
1155  const primitiveMesh& mesh,
1156  const vectorField& faceAreas,
1157  const labelList& checkFaces,
1158  const labelList& affectedCells,
1159  labelHashSet* setPtr
1160 )
1161 {
1162  const cellList& cells = mesh.cells();
1163 
1164  scalar minDet = GREAT;
1165  scalar sumDet = 0.0;
1166  label nSumDet = 0;
1167  label nWarnDet = 0;
1168 
1169  forAll(affectedCells, i)
1170  {
1171  const cell& cFaces = cells[affectedCells[i]];
1172 
1173  tensor areaSum(Zero);
1174  scalar magAreaSum = 0;
1175 
1176  forAll(cFaces, cFacei)
1177  {
1178  label facei = cFaces[cFacei];
1179 
1180  scalar magArea = mag(faceAreas[facei]);
1181 
1182  magAreaSum += magArea;
1183  areaSum += faceAreas[facei]*(faceAreas[facei]/magArea);
1184  }
1185 
1186  scalar scaledDet = det(areaSum/magAreaSum)/0.037037037037037;
1187 
1188  minDet = min(minDet, scaledDet);
1189  sumDet += scaledDet;
1190  nSumDet++;
1191 
1192  if (scaledDet < warnDet)
1193  {
1194  if (setPtr)
1195  {
1196  // Insert all faces of the cell.
1197  forAll(cFaces, cFacei)
1198  {
1199  label facei = cFaces[cFacei];
1200  setPtr->insert(facei);
1201  }
1202  }
1203  nWarnDet++;
1204  }
1205  }
1206 
1207  reduce(minDet, minOp<scalar>());
1208  reduce(sumDet, sumOp<scalar>());
1209  reduce(nSumDet, sumOp<label>());
1210  reduce(nWarnDet, sumOp<label>());
1211 
1212  if (report)
1213  {
1214  if (nSumDet > 0)
1215  {
1216  Info<< "Cell determinant (1 = uniform cube) : average = "
1217  << sumDet / nSumDet << " min = " << minDet << endl;
1218  }
1219 
1220  if (nWarnDet > 0)
1221  {
1222  Info<< "There are " << nWarnDet
1223  << " cells with determinant < " << warnDet << '.' << nl
1224  << endl;
1225  }
1226  else
1227  {
1228  Info<< "All faces have determinant > " << warnDet << '.' << nl
1229  << endl;
1230  }
1231  }
1232 
1233  if (nWarnDet > 0)
1234  {
1235  if (report)
1236  {
1238  << nWarnDet << " cells with determinant < " << warnDet
1239  << " found.\n"
1240  << endl;
1241  }
1242 
1243  return true;
1244  }
1245 
1246  return false;
1247 }
1248 
1249 
1252  const bool report,
1253  const scalar orthWarn,
1254  const labelList& checkFaces,
1255  labelHashSet* setPtr
1256 ) const
1257 {
1258  return checkFaceDotProduct
1259  (
1260  report,
1261  orthWarn,
1262  mesh_,
1263  cellCentres_,
1264  faceAreas_,
1265  checkFaces,
1266  setPtr
1267  );
1268 }
1269 
1270 
1273  const bool report,
1274  const scalar minPyrVol,
1275  const pointField& p,
1276  const labelList& checkFaces,
1277  labelHashSet* setPtr
1278 ) const
1279 {
1280  return checkFacePyramids
1281  (
1282  report,
1283  minPyrVol,
1284  mesh_,
1285  cellCentres_,
1286  p,
1287  checkFaces,
1288  setPtr
1289  );
1290 }
1291 
1292 
1295  const bool report,
1296  const scalar internalSkew,
1297  const scalar boundarySkew,
1298  const labelList& checkFaces,
1299  labelHashSet* setPtr
1300 ) const
1301 {
1302  return checkFaceSkewness
1303  (
1304  report,
1305  internalSkew,
1306  boundarySkew,
1307  mesh_,
1308  cellCentres_,
1309  faceCentres_,
1310  faceAreas_,
1311  checkFaces,
1312  setPtr
1313  );
1314 }
1315 
1316 
1319  const bool report,
1320  const scalar warnWeight,
1321  const labelList& checkFaces,
1322  labelHashSet* setPtr
1323 ) const
1324 {
1325  return checkFaceWeights
1326  (
1327  report,
1328  warnWeight,
1329  mesh_,
1330  cellCentres_,
1331  faceCentres_,
1332  faceAreas_,
1333  checkFaces,
1334  setPtr
1335  );
1336 }
1337 
1338 
1341  const bool report,
1342  const scalar maxDeg,
1343  const pointField& p,
1344  const labelList& checkFaces,
1345  labelHashSet* setPtr
1346 ) const
1347 {
1348  return checkFaceAngles
1349  (
1350  report,
1351  maxDeg,
1352  mesh_,
1353  faceAreas_,
1354  p,
1355  checkFaces,
1356  setPtr
1357  );
1358 }
1359 
1360 
1361 //bool Foam::primitiveMeshGeometry::checkFaceFlatness
1362 //(
1363 // const bool report,
1364 // const scalar warnFlatness,
1365 // const pointField& p,
1366 // const labelList& checkFaces,
1367 // labelHashSet* setPtr
1368 //) const
1369 //{
1370 // return checkFaceFlatness
1371 // (
1372 // report,
1373 // warnFlatness,
1374 // mesh_,
1375 // faceAreas_,
1376 // faceCentres_,
1377 // p,
1378 // checkFaces,
1379 // setPtr
1380 // );
1381 //}
1382 
1383 
1386  const bool report,
1387  const scalar minTwist,
1388  const pointField& p,
1389  const labelList& checkFaces,
1390  labelHashSet* setPtr
1391 ) const
1392 {
1393  return checkFaceTwist
1394  (
1395  report,
1396  minTwist,
1397  mesh_,
1398  faceAreas_,
1399  faceCentres_,
1400  p,
1401  checkFaces,
1402  setPtr
1403  );
1404 }
1405 
1406 
1409  const bool report,
1410  const scalar minArea,
1411  const labelList& checkFaces,
1412  labelHashSet* setPtr
1413 ) const
1414 {
1415  return checkFaceArea
1416  (
1417  report,
1418  minArea,
1419  mesh_,
1420  faceAreas_,
1421  checkFaces,
1422  setPtr
1423  );
1424 }
1425 
1426 
1429  const bool report,
1430  const scalar warnDet,
1431  const labelList& checkFaces,
1432  const labelList& affectedCells,
1433  labelHashSet* setPtr
1434 ) const
1435 {
1436  return checkCellDeterminant
1437  (
1438  report,
1439  warnDet,
1440  mesh_,
1441  faceAreas_,
1442  checkFaces,
1443  affectedCells,
1444  setPtr
1445  );
1446 }
1447 
1448 
1449 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::Tensor< scalar >
Foam::maxOp
Definition: ops.H:223
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::primitiveMeshGeometry::affectedCells
labelList affectedCells(const labelList &changedFaces) const
Helper function: get affected cells from faces.
Definition: primitiveMeshGeometry.C:192
s
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))
Definition: gmvOutputSpray.H:25
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::primitiveMesh::faces
virtual const faceList & faces() const =0
Return faces.
Foam::minOp
Definition: ops.H:224
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::primitiveMeshGeometry::checkFacePyramids
static bool checkFacePyramids(const bool report, const scalar minPyrVol, const primitiveMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, labelHashSet *)
Definition: primitiveMeshGeometry.C:400
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
triPointRef.H
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::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMeshGeometry::checkCellDeterminant
static bool checkCellDeterminant(const bool report, const scalar minDet, const primitiveMesh &, const vectorField &faceAreas, const labelList &checkFaces, const labelList &affectedCells, labelHashSet *setPtr)
Definition: primitiveMeshGeometry.C:1152
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
correct
fvOptions correct(rho)
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
SeriousErrorInFunction
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
Definition: messageStream.H:306
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::primitiveMeshGeometry::checkFaceAngles
static bool checkFaceAngles(const bool report, const scalar maxDeg, const primitiveMesh &mesh, const vectorField &faceAreas, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Definition: primitiveMeshGeometry.C:731
Foam::primitiveMeshGeometry::checkFaceSkewness
static bool checkFaceSkewness(const bool report, const scalar internalSkew, const scalar boundarySkew, const primitiveMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Definition: primitiveMeshGeometry.C:510
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::primitiveMeshGeometry::correct
void correct()
Take over properties from mesh.
Definition: primitiveMeshGeometry.C:236
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:103
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::primitiveMeshGeometry::primitiveMeshGeometry
primitiveMeshGeometry(const primitiveMesh &)
Construct from mesh.
Definition: primitiveMeshGeometry.C:221
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::primitiveMeshGeometry::checkFaceArea
static bool checkFaceArea(const bool report, const scalar minArea, const primitiveMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Definition: primitiveMeshGeometry.C:1093
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMeshGeometry::checkFaceDotProduct
static bool checkFaceDotProduct(const bool report, const scalar orthWarn, const primitiveMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Definition: primitiveMeshGeometry.C:259
Foam::primitiveMeshGeometry::checkFaceWeights
static bool checkFaceWeights(const bool report, const scalar warnWeight, const primitiveMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Definition: primitiveMeshGeometry.C:647
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
primitiveMeshGeometry.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::det
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:62
Foam::primitiveMeshGeometry::checkFaceTwist
static bool checkFaceTwist(const bool report, const scalar minTwist, const primitiveMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Definition: primitiveMeshGeometry.C:989
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::pyramidPointFaceRef
pyramid< point, const point &, const face & > pyramidPointFaceRef
A pyramid using referred points and faces.
Definition: pyramidPointFaceRef.H:48
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
pyramidPointFaceRef.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78