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