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