polyMeshCheck.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 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 "polyMesh.H"
30 #include "polyMeshTools.H"
31 #include "unitConversion.H"
32 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 bool Foam::polyMesh::checkFaceOrthogonality
37 (
38  const vectorField& fAreas,
39  const vectorField& cellCtrs,
40  const bool report,
41  const bool detailedReport,
42  labelHashSet* setPtr
43 ) const
44 {
45  DebugInFunction << "Checking mesh non-orthogonality" << endl;
46 
47  const labelList& own = faceOwner();
48  const labelList& nei = faceNeighbour();
49 
50  // Calculate orthogonality for all internal and coupled boundary faces
51  // (1 for uncoupled boundary faces)
52  tmp<scalarField> tortho = polyMeshTools::faceOrthogonality
53  (
54  *this,
55  fAreas,
56  cellCtrs
57  );
58  const scalarField& ortho = tortho.ref();
59 
60  // Severe nonorthogonality threshold
61  const scalar severeNonorthogonalityThreshold =
63 
64 
65  scalar minDDotS = GREAT;
66  scalar sumDDotS = 0.0;
67  label nSummed = 0;
68  label severeNonOrth = 0;
69  label errorNonOrth = 0;
70 
71 
72  // Statistics only for internal and masters of coupled faces
73  bitSet isMasterFace(syncTools::getInternalOrMasterFaces(*this));
74 
75  forAll(ortho, facei)
76  {
77  if (ortho[facei] < severeNonorthogonalityThreshold)
78  {
79  if (ortho[facei] > SMALL)
80  {
81  if (setPtr)
82  {
83  setPtr->insert(facei);
84  }
85 
86  severeNonOrth++;
87  }
88  else
89  {
90  // Error : non-ortho too large
91  if (setPtr)
92  {
93  setPtr->insert(facei);
94  }
95  if (detailedReport && errorNonOrth == 0)
96  {
97  // Non-orthogonality greater than 90 deg
99  << "Severe non-orthogonality for face "
100  << facei
101  << " between cells " << own[facei]
102  << " and " << nei[facei]
103  << ": Angle = "
104  << radToDeg(::acos(min(1.0, max(-1.0, ortho[facei]))))
105  << " deg." << endl;
106  }
107 
108  errorNonOrth++;
109  }
110  }
111 
112  if (isMasterFace.test(facei))
113  {
114  minDDotS = min(minDDotS, ortho[facei]);
115  sumDDotS += ortho[facei];
116  nSummed++;
117  }
118  }
119 
120  reduce(minDDotS, minOp<scalar>());
121  reduce(sumDDotS, sumOp<scalar>());
122  reduce(nSummed, sumOp<label>());
123  reduce(severeNonOrth, sumOp<label>());
124  reduce(errorNonOrth, sumOp<label>());
125 
126  if (debug || report)
127  {
128  if (nSummed > 0)
129  {
130  if (debug || report)
131  {
132  Info<< " Mesh non-orthogonality Max: "
133  << radToDeg(::acos(min(1.0, max(-1.0, minDDotS))))
134  << " average: "
135  << radToDeg(::acos(min(1.0, max(-1.0, sumDDotS/nSummed))))
136  << endl;
137  }
138  }
139 
140  if (severeNonOrth > 0)
141  {
142  Info<< " *Number of severely non-orthogonal (> "
143  << primitiveMesh::nonOrthThreshold_ << " degrees) faces: "
144  << severeNonOrth << "." << endl;
145  }
146  }
147 
148  if (errorNonOrth > 0)
149  {
150  if (debug || report)
151  {
152  Info<< " ***Number of non-orthogonality errors: "
153  << errorNonOrth << "." << endl;
154  }
155 
156  return true;
157  }
158 
159  if (debug || report)
160  {
161  Info<< " Non-orthogonality check OK." << endl;
162  }
163 
164  return false;
165 }
166 
167 
168 bool Foam::polyMesh::checkFaceSkewness
169 (
170  const pointField& points,
171  const vectorField& fCtrs,
172  const vectorField& fAreas,
173  const vectorField& cellCtrs,
174  const bool report,
175  const bool detailedReport,
176  labelHashSet* setPtr
177 ) const
178 {
179  DebugInFunction << "Checking face skewness" << endl;
180 
181  const labelList& own = faceOwner();
182  const labelList& nei = faceNeighbour();
183 
184  // Warn if the skew correction vector is more than skewWarning times
185  // larger than the face area vector
186 
187  tmp<scalarField> tskew = polyMeshTools::faceSkewness
188  (
189  *this,
190  points,
191  fCtrs,
192  fAreas,
193  cellCtrs
194  );
195  const scalarField& skew = tskew.ref();
196 
197  scalar maxSkew = max(skew);
198  label nWarnSkew = 0;
199 
200  // Statistics only for all faces except slave coupled faces
201  bitSet isMasterFace(syncTools::getMasterFaces(*this));
202 
203  forAll(skew, facei)
204  {
205  // Check if the skewness vector is greater than the PN vector.
206  // This does not cause trouble but is a good indication of a poor mesh.
207  if (skew[facei] > skewThreshold_)
208  {
209  if (setPtr)
210  {
211  setPtr->insert(facei);
212  }
213  if (detailedReport && nWarnSkew == 0)
214  {
215  // Non-orthogonality greater than 90 deg
216  if (isInternalFace(facei))
217  {
219  << "Severe skewness " << skew[facei]
220  << " for face " << facei
221  << " between cells " << own[facei]
222  << " and " << nei[facei];
223  }
224  else
225  {
227  << "Severe skewness " << skew[facei]
228  << " for boundary face " << facei
229  << " on cell " << own[facei];
230  }
231  }
232 
233  if (isMasterFace.test(facei))
234  {
235  ++nWarnSkew;
236  }
237  }
238  }
239 
240  reduce(maxSkew, maxOp<scalar>());
241  reduce(nWarnSkew, sumOp<label>());
242 
243  if (nWarnSkew > 0)
244  {
245  if (debug || report)
246  {
247  Info<< " ***Max skewness = " << maxSkew
248  << ", " << nWarnSkew << " highly skew faces detected"
249  " which may impair the quality of the results"
250  << endl;
251  }
252 
253  return true;
254  }
255 
256  if (debug || report)
257  {
258  Info<< " Max skewness = " << maxSkew << " OK." << endl;
259  }
260 
261  return false;
262 }
263 
264 
265 bool Foam::polyMesh::checkEdgeAlignment
266 (
267  const pointField& p,
268  const bool report,
269  const Vector<label>& directions,
270  labelHashSet* setPtr
271 ) const
272 {
273  // Check 1D/2Dness of edges. Gets passed the non-empty directions and
274  // checks all edges in the mesh whether they:
275  // - have no component in a non-empty direction or
276  // - are only in a single non-empty direction.
277  // Empty direction info is passed in as a vector of labels (synchronised)
278  // which are 1 if the direction is non-empty, 0 if it is.
279 
280  DebugInFunction << "Checking edge alignment" << endl;
281 
282  label nDirs = 0;
283  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
284  {
285  if (directions[cmpt] == 1)
286  {
287  nDirs++;
288  }
289  else if (directions[cmpt] != 0)
290  {
292  << "directions should contain 0 or 1 but is now " << directions
293  << exit(FatalError);
294  }
295  }
296 
297  if (nDirs == vector::nComponents)
298  {
299  return false;
300  }
301 
302 
303  const faceList& fcs = faces();
304 
305  EdgeMap<label> edgesInError;
306 
307  forAll(fcs, facei)
308  {
309  const face& f = fcs[facei];
310 
311  forAll(f, fp)
312  {
313  label p0 = f[fp];
314  label p1 = f.nextLabel(fp);
315  if (p0 < p1)
316  {
317  vector d(p[p1]-p[p0]);
318  scalar magD = mag(d);
319 
320  if (magD > ROOTVSMALL)
321  {
322  d /= magD;
323 
324  // Check how many empty directions are used by the edge.
325  label nEmptyDirs = 0;
326  label nNonEmptyDirs = 0;
327  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
328  {
329  if (mag(d[cmpt]) > 1e-6)
330  {
331  if (directions[cmpt] == 0)
332  {
333  nEmptyDirs++;
334  }
335  else
336  {
337  nNonEmptyDirs++;
338  }
339  }
340  }
341 
342  if (nEmptyDirs == 0)
343  {
344  // Purely in ok directions.
345  }
346  else if (nEmptyDirs == 1)
347  {
348  // Ok if purely in empty directions.
349  if (nNonEmptyDirs > 0)
350  {
351  edgesInError.insert(edge(p0, p1), facei);
352  }
353  }
354  else if (nEmptyDirs > 1)
355  {
356  // Always an error
357  edgesInError.insert(edge(p0, p1), facei);
358  }
359  }
360  }
361  }
362  }
363 
364  label nErrorEdges = returnReduce(edgesInError.size(), sumOp<label>());
365 
366  if (nErrorEdges > 0)
367  {
368  if (debug || report)
369  {
370  Info<< " ***Number of edges not aligned with or perpendicular to "
371  << "non-empty directions: " << nErrorEdges << endl;
372  }
373 
374  if (setPtr)
375  {
376  setPtr->resize(2*edgesInError.size());
377  forAllConstIters(edgesInError, iter)
378  {
379  setPtr->insert(iter.key()[0]);
380  setPtr->insert(iter.key()[1]);
381  }
382  }
383 
384  return true;
385  }
386 
387  if (debug || report)
388  {
389  Info<< " All edges aligned with or perpendicular to "
390  << "non-empty directions." << endl;
391  }
392 
393  return false;
394 }
395 
396 
397 bool Foam::polyMesh::checkCellDeterminant
398 (
399  const vectorField& faceAreas,
400  const bool report,
401  labelHashSet* setPtr,
402  const Vector<label>& meshD
403 ) const
404 {
405  const scalar warnDet = 1e-3;
406 
407  DebugInFunction << "Checking for under-determined cells" << endl;
408 
409  tmp<scalarField> tcellDeterminant = primitiveMeshTools::cellDeterminant
410  (
411  *this,
412  meshD,
413  faceAreas,
415  );
416  scalarField& cellDeterminant = tcellDeterminant.ref();
417 
418 
419  label nErrorCells = 0;
420  scalar minDet = min(cellDeterminant);
421  scalar sumDet = sum(cellDeterminant);
422 
423  forAll(cellDeterminant, celli)
424  {
425  if (cellDeterminant[celli] < warnDet)
426  {
427  if (setPtr)
428  {
429  setPtr->insert(celli);
430  }
431 
432  nErrorCells++;
433  }
434  }
435 
436  reduce(nErrorCells, sumOp<label>());
437  reduce(minDet, minOp<scalar>());
438  reduce(sumDet, sumOp<scalar>());
439  label nSummed = returnReduce(cellDeterminant.size(), sumOp<label>());
440 
441  if (debug || report)
442  {
443  if (nSummed > 0)
444  {
445  Info<< " Cell determinant (wellposedness) : minimum: " << minDet
446  << " average: " << sumDet/nSummed
447  << endl;
448  }
449  }
450 
451  if (nErrorCells > 0)
452  {
453  if (debug || report)
454  {
455  Info<< " ***Cells with small determinant (< "
456  << warnDet << ") found, number of cells: "
457  << nErrorCells << endl;
458  }
459 
460  return true;
461  }
462 
463  if (debug || report)
464  {
465  Info<< " Cell determinant check OK." << endl;
466  }
467 
468  return false;
469 }
470 
471 
472 bool Foam::polyMesh::checkFaceWeight
473 (
474  const vectorField& fCtrs,
475  const vectorField& fAreas,
476  const vectorField& cellCtrs,
477  const bool report,
478  const scalar minWeight,
479  labelHashSet* setPtr
480 ) const
481 {
482  DebugInFunction << "Checking for low face interpolation weights" << endl;
483 
484  tmp<scalarField> tfaceWght = polyMeshTools::faceWeights
485  (
486  *this,
487  fCtrs,
488  fAreas,
489  cellCtrs
490  );
491  scalarField& faceWght = tfaceWght.ref();
492 
493 
494  label nErrorFaces = 0;
495  scalar minDet = GREAT;
496  scalar sumDet = 0.0;
497  label nSummed = 0;
498 
499  // Statistics only for internal and masters of coupled faces
500  bitSet isMasterFace(syncTools::getInternalOrMasterFaces(*this));
501 
502  forAll(faceWght, facei)
503  {
504  if (faceWght[facei] < minWeight)
505  {
506  // Note: insert both sides of coupled faces
507  if (setPtr)
508  {
509  setPtr->insert(facei);
510  }
511 
512  nErrorFaces++;
513  }
514 
515  // Note: statistics only on master of coupled faces
516  if (isMasterFace.test(facei))
517  {
518  minDet = min(minDet, faceWght[facei]);
519  sumDet += faceWght[facei];
520  nSummed++;
521  }
522  }
523 
524  reduce(nErrorFaces, sumOp<label>());
525  reduce(minDet, minOp<scalar>());
526  reduce(sumDet, sumOp<scalar>());
527  reduce(nSummed, sumOp<label>());
528 
529  if (debug || report)
530  {
531  if (nSummed > 0)
532  {
533  Info<< " Face interpolation weight : minimum: " << minDet
534  << " average: " << sumDet/nSummed
535  << endl;
536  }
537  }
538 
539  if (nErrorFaces > 0)
540  {
541  if (debug || report)
542  {
543  Info<< " ***Faces with small interpolation weight (< " << minWeight
544  << ") found, number of faces: "
545  << nErrorFaces << endl;
546  }
547 
548  return true;
549  }
550 
551  if (debug || report)
552  {
553  Info<< " Face interpolation weight check OK." << endl;
554  }
555 
556  return false;
557 }
558 
559 
560 bool Foam::polyMesh::checkVolRatio
561 (
562  const scalarField& cellVols,
563  const bool report,
564  const scalar minRatio,
565  labelHashSet* setPtr
566 ) const
567 {
568  DebugInFunction << "Checking for volume ratio < " << minRatio << endl;
569 
570  tmp<scalarField> tvolRatio = polyMeshTools::volRatio(*this, cellVols);
571  scalarField& volRatio = tvolRatio.ref();
572 
573 
574  label nErrorFaces = 0;
575  scalar minDet = GREAT;
576  scalar sumDet = 0.0;
577  label nSummed = 0;
578 
579  // Statistics only for internal and masters of coupled faces
580  bitSet isMasterFace(syncTools::getInternalOrMasterFaces(*this));
581 
582  forAll(volRatio, facei)
583  {
584  if (volRatio[facei] < minRatio)
585  {
586  // Note: insert both sides of coupled faces
587  if (setPtr)
588  {
589  setPtr->insert(facei);
590  }
591 
592  nErrorFaces++;
593  }
594 
595  // Note: statistics only on master of coupled faces
596  if (isMasterFace.test(facei))
597  {
598  minDet = min(minDet, volRatio[facei]);
599  sumDet += volRatio[facei];
600  nSummed++;
601  }
602  }
603 
604  reduce(nErrorFaces, sumOp<label>());
605  reduce(minDet, minOp<scalar>());
606  reduce(sumDet, sumOp<scalar>());
607  reduce(nSummed, sumOp<label>());
608 
609  if (debug || report)
610  {
611  if (nSummed > 0)
612  {
613  Info<< " Face volume ratio : minimum: " << minDet
614  << " average: " << sumDet/nSummed
615  << endl;
616  }
617  }
618 
619  if (nErrorFaces > 0)
620  {
621  if (debug || report)
622  {
623  Info<< " ***Faces with small volume ratio (< " << minRatio
624  << ") found, number of faces: "
625  << nErrorFaces << endl;
626  }
627 
628  return true;
629  }
630 
631  if (debug || report)
632  {
633  Info<< " Face volume ratio check OK." << endl;
634  }
635 
636  return false;
637 }
638 
639 
640 bool Foam::polyMesh::checkFaceOrthogonality
641 (
642  const bool report,
643  labelHashSet* setPtr
644 ) const
645 {
646  return checkFaceOrthogonality
647  (
648  faceAreas(),
649  cellCentres(),
650  report,
651  false, // detailedReport
652  setPtr
653  );
654 }
655 
656 
657 bool Foam::polyMesh::checkFaceSkewness
658 (
659  const bool report,
660  labelHashSet* setPtr
661 ) const
662 {
663  return checkFaceSkewness
664  (
665  points(),
666  faceCentres(),
667  faceAreas(),
668  cellCentres(),
669  report,
670  false, // detailedReport
671  setPtr
672  );
673 }
674 
675 
676 bool Foam::polyMesh::checkEdgeAlignment
677 (
678  const bool report,
679  const Vector<label>& directions,
680  labelHashSet* setPtr
681 ) const
682 {
683  return checkEdgeAlignment
684  (
685  points(),
686  report,
687  directions,
688  setPtr
689  );
690 }
691 
692 
693 bool Foam::polyMesh::checkCellDeterminant
694 (
695  const bool report,
696  labelHashSet* setPtr
697 ) const
698 {
699  return checkCellDeterminant
700  (
701  faceAreas(),
702  report,
703  setPtr,
704  geometricD()
705  );
706 }
707 
708 
709 bool Foam::polyMesh::checkFaceWeight
710 (
711  const bool report,
712  const scalar minWeight,
713  labelHashSet* setPtr
714 ) const
715 {
716  return checkFaceWeight
717  (
718  faceCentres(),
719  faceAreas(),
720  cellCentres(),
721  report,
722  minWeight,
723  setPtr
724  );
725 }
726 
727 
728 bool Foam::polyMesh::checkVolRatio
729 (
730  const bool report,
731  const scalar minRatio,
732  labelHashSet* setPtr
733 ) const
734 {
735  return checkVolRatio(cellVolumes(), report, minRatio, setPtr);
736 }
737 
738 
740 (
741  const pointField& newPoints,
742  const bool report,
743  const bool detailedReport
744 ) const
745 {
746  if (debug || report)
747  {
748  Pout<< "bool polyMesh::checkMeshMotion("
749  << "const pointField&, const bool, const bool) const: "
750  << "checking mesh motion" << endl;
751  }
752 
753  vectorField fCtrs(nFaces());
754  vectorField fAreas(nFaces());
755 
757  (
758  *this,
759  newPoints,
760  fCtrs,
761  fAreas
762  );
763 
764  // Check cell volumes and calculate new cell centres
765  vectorField cellCtrs(nCells());
766  scalarField cellVols(nCells());
767 
769  (
770  *this,
771  fCtrs,
772  fAreas,
773  cellCtrs,
774  cellVols
775  );
776 
777  // Check cell volumes
778  bool error = checkCellVolumes
779  (
780  cellVols, // vols
781  report, // report
782  detailedReport, // detailedReport
783  nullptr // setPtr
784  );
785 
786 
787  // Check face areas
788  bool areaError = checkFaceAreas
789  (
790  fAreas,
791  report, // report
792  detailedReport, // detailedReport,
793  nullptr // setPtr
794  );
795  error = error || areaError;
796 
797 
798  // Check pyramid volumes
799  bool pyrVolError = checkFacePyramids
800  (
801  newPoints,
802  cellCtrs,
803  report, // report,
804  detailedReport, // detailedReport,
805  -SMALL, // minPyrVol
806  nullptr // setPtr
807  );
808  error = error || pyrVolError;
809 
810 
811  // Check face non-orthogonality
812  bool nonOrthoError = checkFaceOrthogonality
813  (
814  fAreas,
815  cellCtrs,
816  report, // report
817  detailedReport, // detailedReport
818  nullptr // setPtr
819  );
820  error = error || nonOrthoError;
821 
822 
823  if (!error && (debug || report))
824  {
825  Pout<< "Mesh motion check OK." << endl;
826  }
827 
828  return error;
829 }
830 
831 
832 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::primitiveMeshTools::cellDeterminant
static tmp< scalarField > cellDeterminant(const primitiveMesh &mesh, const Vector< label > &directions, const vectorField &faceAreas, const bitSet &internalOrCoupledFace)
Generate cell determinant field. Normalised to 1 for an internal cube.
Definition: primitiveMeshTools.C:628
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::syncTools::getInternalOrCoupledFaces
static bitSet getInternalOrCoupledFaces(const polyMesh &mesh)
Get per face whether it is internal or coupled.
Definition: syncTools.C:176
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:138
Foam::polyMeshTools::volRatio
static tmp< scalarField > volRatio(const polyMesh &mesh, const scalarField &vol)
Generate volume ratio field.
Definition: polyMeshTools.C:235
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::syncTools::getMasterFaces
static bitSet getMasterFaces(const polyMesh &mesh)
Definition: syncTools.C:126
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::polyMeshTools::faceSkewness
static tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
Definition: polyMeshTools.C:92
polyMesh.H
Foam::directions
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:67
Foam::HashSet< label, Hash< label > >
syncTools.H
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
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
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::primitiveMeshTools::makeCellCentresAndVols
static void makeCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols)
Calculate cell centres and volumes from face properties.
Definition: primitiveMeshTools.C:107
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
cellVols
const scalarField & cellVols
Definition: temperatureAndPressureVariables.H:51
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::FatalError
error FatalError
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::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
polyMeshTools.H
Foam::primitiveMeshTools::makeFaceCentresAndAreas
static void makeFaceCentresAndAreas(const primitiveMesh &mesh, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Calculate face centres and areas.
Definition: primitiveMeshTools.C:38
Foam::primitiveMesh::nonOrthThreshold_
static scalar nonOrthThreshold_
Non-orthogonality warning threshold in deg.
Definition: primitiveMesh.H:244
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::Vector< label >
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::polyMesh::checkMeshMotion
virtual bool checkMeshMotion(const pointField &newPoints, const bool report=false, const bool detailedReport=false) const
Check mesh motion for correctness given motion points.
Definition: polyMeshCheck.C:740
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::polyMeshTools::faceWeights
static tmp< scalarField > faceWeights(const polyMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate interpolation factors field.
Definition: polyMeshTools.C:177
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::error
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:73
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::polyMeshTools::faceOrthogonality
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for.
Definition: polyMeshTools.C:37
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::syncTools::getInternalOrMasterFaces
static bitSet getInternalOrMasterFaces(const polyMesh &mesh)
Definition: syncTools.C:148