52 const bitSet& internalOrCoupledFaces
61 solveScalar sumMagClosedBoundary(0);
65 if (!internalOrCoupledFaces.
size() || !internalOrCoupledFaces[facei])
67 sumClosed += areas[facei];
68 sumMagClosedBoundary +=
mag(areas[facei]);
75 solveVector openness = sumClosed/(sumMagClosedBoundary + VSMALL);
81 Info<<
" ***Boundary openness " << openness
82 <<
" possible hole in boundary description."
91 Info<<
" Boundary openness " << openness <<
" OK."
106 const Vector<label>& meshD
114 label nErrorClosed = 0;
118 const cell& curCell = c[cI];
120 if (
min(curCell) < 0 ||
max(curCell) > nFaces())
131 if (nErrorClosed > 0)
135 Info<<
" ***Cells with invalid face labels found, number of cells "
136 << nErrorClosed <<
endl;
156 scalar maxOpennessCell =
max(openness);
158 scalar maxAspectRatio =
max(aspectRatio);
163 if (openness[celli] > closedThreshold_)
167 setPtr->insert(celli);
173 if (aspectRatio[celli] > aspectThreshold_)
177 aspectSetPtr->insert(celli);
195 Info<<
" ***Open cells found, max cell openness: "
196 << maxOpennessCell <<
", number of open cells " << nOpen
207 Info<<
" ***High aspect ratio cells found, Max aspect ratio: "
209 <<
", number of cells " << nAspect
218 Info<<
" Max cell openness = " << maxOpennessCell <<
" OK." <<
nl
219 <<
" Max aspect ratio = " << maxAspectRatio <<
" OK."
231 const bool detailedReport,
239 scalar minArea = GREAT;
240 scalar maxArea = -GREAT;
242 forAll(magFaceAreas, facei)
244 if (magFaceAreas[facei] < VSMALL)
248 setPtr->insert(facei);
252 if (isInternalFace(facei))
254 Pout<<
"Zero or negative face area detected for "
255 <<
"internal face "<< facei <<
" between cells "
256 << faceOwner()[facei] <<
" and "
257 << faceNeighbour()[facei]
258 <<
". Face area magnitude = " << magFaceAreas[facei]
263 Pout<<
"Zero or negative face area detected for "
264 <<
"boundary face " << facei <<
" next to cell "
265 << faceOwner()[facei] <<
". Face area magnitude = "
266 << magFaceAreas[facei] <<
endl;
271 minArea =
min(minArea, magFaceAreas[facei]);
272 maxArea =
max(maxArea, magFaceAreas[facei]);
278 if (minArea < VSMALL)
282 Info<<
" ***Zero or negative face area detected. "
283 "Minimum area: " << minArea <<
endl;
291 Info<<
" Minimum face area = " << minArea
292 <<
". Maximum face area = " << maxArea
293 <<
". Face area magnitudes OK." <<
endl;
304 const bool detailedReport,
310 scalar minVolume = GREAT;
311 scalar maxVolume = -GREAT;
313 label nNegVolCells = 0;
317 if (vols[celli] < VSMALL)
321 setPtr->insert(celli);
325 Pout<<
"Zero or negative cell volume detected for cell "
326 << celli <<
". Volume = " << vols[celli] <<
endl;
332 minVolume =
min(minVolume, vols[celli]);
333 maxVolume =
max(maxVolume, vols[celli]);
340 if (minVolume < VSMALL)
344 Info<<
" ***Zero or negative cell volume detected. "
345 <<
"Minimum negative volume: " << minVolume
346 <<
", Number of negative volume cells: " << nNegVolCells
355 Info<<
" Min volume = " << minVolume
356 <<
". Max volume = " << maxVolume
357 <<
". Total volume = " <<
gSum(vols)
358 <<
". Cell volumes OK." <<
endl;
384 const scalar severeNonorthogonalityThreshold =
387 scalar minDDotS =
min(ortho);
389 scalar sumDDotS =
sum(ortho);
391 label severeNonOrth = 0;
393 label errorNonOrth = 0;
398 if (ortho[facei] < severeNonorthogonalityThreshold)
400 if (ortho[facei] > SMALL)
404 setPtr->insert(facei);
413 setPtr->insert(facei);
434 Info<<
" Mesh non-orthogonality Max: "
441 if (severeNonOrth > 0)
443 Info<<
" *Number of severely non-orthogonal faces: "
444 << severeNonOrth <<
"." <<
endl;
448 if (errorNonOrth > 0)
452 Info<<
" ***Number of non-orthogonality errors: "
453 << errorNonOrth <<
"." <<
endl;
473 const bool detailedReport,
474 const scalar minPyrVol,
497 label nErrorPyrs = 0;
501 if (ownPyrVol[facei] < minPyrVol)
505 setPtr->insert(facei);
509 Pout<<
"Negative pyramid volume: " << ownPyrVol[facei]
510 <<
" for face " << facei <<
" " <<
f[facei]
511 <<
" and owner cell: " << own[facei] <<
endl
512 <<
"Owner cell vertex labels: "
513 <<
cells()[own[facei]].labels(faces())
520 if (isInternalFace(facei))
522 if (neiPyrVol[facei] < minPyrVol)
526 setPtr->insert(facei);
530 Pout<<
"Negative pyramid volume: " << neiPyrVol[facei]
531 <<
" for face " << facei <<
" " <<
f[facei]
532 <<
" and neighbour cell: " << nei[facei] <<
nl
533 <<
"Neighbour cell vertex labels: "
534 <<
cells()[nei[facei]].labels(faces())
548 Info<<
" ***Error in face pyramids: "
549 << nErrorPyrs <<
" faces are incorrectly oriented."
590 scalar maxSkew =
max(skewness);
597 if (skewness[facei] > skewThreshold_)
601 setPtr->insert(facei);
615 Info<<
" ***Max skewness = " << maxSkew
616 <<
", " << nWarnSkew <<
" highly skew faces detected"
617 " which may impair the quality of the results"
626 Info<<
" Max skewness = " << maxSkew <<
" OK." <<
endl;
644 if (maxDeg < -SMALL || maxDeg > 180+SMALL)
647 <<
"maxDeg should be [0..180] but is now " << maxDeg
663 scalar maxEdgeSin =
max(faceAngles);
669 if (faceAngles[facei] > SMALL)
675 setPtr->insert(facei);
685 scalar maxConcaveDegr =
690 Info<<
" *There are " << nConcave
691 <<
" faces with concave angles between consecutive"
692 <<
" edges. Max concave angle = " << maxConcaveDegr
693 <<
" degrees." <<
endl;
714 const scalar warnFlatness,
720 if (warnFlatness < 0 || warnFlatness > 1)
723 <<
"warnFlatness should be [0..1] but is " << warnFlatness
740 scalar minFlatness = GREAT;
741 scalar sumFlatness = 0;
745 forAll(faceFlatness, facei)
747 if (fcs[facei].size() > 3 && magAreas[facei] > VSMALL)
749 sumFlatness += faceFlatness[facei];
752 minFlatness =
min(minFlatness, faceFlatness[facei]);
754 if (faceFlatness[facei] < warnFlatness)
760 setPtr->insert(facei);
777 Info<<
" Face flatness (1 = flat, 0 = butterfly) : min = "
778 << minFlatness <<
" average = " << sumFlatness / nSummed
788 Info<<
" *There are " << nWarped
789 <<
" faces with ratio between projected and actual area < "
790 << warnFlatness <<
endl;
792 Info<<
" Minimum ratio (minimum flatness, maximum warpage) = "
793 << minFlatness <<
endl;
821 label nConcaveCells = 0;
825 const cell& cFaces = c[celli];
827 bool concave =
false;
836 label fI = cFaces[i];
838 const point& fC = fCentres[fI];
842 fN /=
max(
mag(fN), VSMALL);
846 if (fOwner[fI] != celli)
858 label fJ = cFaces[j];
860 const point& pt = fCentres[fJ];
870 pC /=
max(
mag(pC), VSMALL);
872 if ((pC & fN) > -planarCosAngle_)
880 setPtr->insert(celli);
894 if (nConcaveCells > 0)
898 Info<<
" ***Concave cells (using face planes) found,"
899 <<
" number of cells: " << nConcaveCells <<
endl;
928 label internal = nInternalFaces();
931 bool hasError =
false;
933 label nMultipleCells = 0;
937 for (label facei = 0; facei < internal; facei++)
939 if (own[facei] >= nei[facei])
945 setPtr->insert(facei);
963 label facei = curFaces[i];
965 if (facei >= nInternalFaces())
972 label nbrCelli = nei[facei];
974 if (nbrCelli == celli)
976 nbrCelli = own[facei];
979 if (celli < nbrCelli)
999 label prevCell = nbr[0];
1000 label prevFace = curFaces[nbr.indices()[0]];
1002 bool hasMultipleFaces =
false;
1004 for (label i = 1; i < nbr.size(); i++)
1006 label thisCell = nbr[i];
1007 label thisFace = curFaces[nbr.indices()[i]];
1014 if (thisCell == prevCell)
1016 hasMultipleFaces =
true;
1020 setPtr->insert(prevFace);
1021 setPtr->insert(thisFace);
1024 else if (thisFace < prevFace)
1030 setPtr->insert(thisFace);
1034 prevCell = thisCell;
1035 prevFace = thisFace;
1038 if (hasMultipleFaces)
1047 if ((
debug || report) && nMultipleCells > 0)
1049 Info<<
" <<Found " << nMultipleCells
1050 <<
" neighbouring cells with multiple inbetween faces." <<
endl;
1055 if (
debug || report)
1057 Info<<
" ***Faces not in upper triangular order." <<
endl;
1063 if (
debug || report)
1080 label nOpenCells = 0;
1089 const edgeList cellEdges = c[celli].edges(
f);
1095 edgeList curFaceEdges =
f[curFaces[facei]].edges();
1097 forAll(curFaceEdges, faceEdgeI)
1099 const edge& curEdge = curFaceEdges[faceEdgeI];
1101 forAll(cellEdges, cellEdgeI)
1103 if (cellEdges[cellEdgeI] == curEdge)
1105 edgeUsage[cellEdgeI]++;
1112 edgeList singleEdges(cellEdges.size());
1113 label nSingleEdges = 0;
1117 if (edgeUsage[edgeI] == 1)
1119 singleEdges[nSingleEdges] = cellEdges[edgeI];
1122 else if (edgeUsage[edgeI] != 2)
1126 setPtr->insert(celli);
1131 if (nSingleEdges > 0)
1135 setPtr->insert(celli);
1146 if (
debug || report)
1148 Info<<
" ***Open cells found, number of cells: " << nOpenCells
1149 <<
". This problem may be fixable using the zipUpMesh utility."
1156 if (
debug || report)
1176 label nErrorFaces = 0;
1180 const face& curFace =
f[fI];
1197 bool inserted = facePoints.insert(curFace[fp]);
1213 if (nErrorFaces > 0)
1215 if (
debug || report)
1217 Info<<
" Faces with invalid vertex labels found, "
1218 <<
" number of faces: " << nErrorFaces <<
endl;
1224 if (
debug || report)
1241 label nFaceErrors = 0;
1242 label nCellErrors = 0;
1248 if (pf[pointi].empty())
1252 setPtr->insert(pointi);
1268 setPtr->insert(pointi);
1278 if (nFaceErrors > 0 || nCellErrors > 0)
1280 if (
debug || report)
1282 Info<<
" ***Unused points found in the mesh, "
1283 "number unused by faces: " << nFaceErrors
1284 <<
" number unused by cells: " << nCellErrors
1291 if (
debug || report)
1303 const Map<label>& nCommonPoints,
1304 label& nBaffleFaces,
1312 label nbFacei = iter.key();
1313 label nCommon = iter();
1315 const face& curFace = faces()[facei];
1316 const face& nbFace = faces()[nbFacei];
1318 if (nCommon == nbFace.size() || nCommon == curFace.size())
1320 if (nbFace.size() != curFace.size())
1331 setPtr->insert(facei);
1332 setPtr->insert(nbFacei);
1352 label nbFacei = iter.key();
1353 label nCommon = iter();
1355 const face& curFace = faces()[facei];
1356 const face& nbFace = faces()[nbFacei];
1361 && nCommon != nbFace.size()
1362 && nCommon != curFace.size()
1368 label nb = nbFace.find(curFace[fp]);
1382 label fpPlus1 = curFace.fcIndex(fp);
1383 label fpMin1 = curFace.rcIndex(fp);
1386 label nbPlus1 = nbFace.fcIndex(nb);
1387 label nbMin1 = nbFace.rcIndex(nb);
1394 if (nbFace[nbPlus1] == curFace[fpPlus1])
1399 else if (nbFace[nbPlus1] == curFace[fpMin1])
1404 else if (nbFace[nbMin1] == curFace[fpMin1])
1424 if (curFp >= curFace.size())
1430 curFp = curFace.size()-1;
1435 if (curNb >= nbFace.size())
1441 curNb = nbFace.size()-1;
1443 }
while (curFace[curFp] == nbFace[curNb]);
1452 for (label commonI = 0; commonI < nCommon; commonI++)
1456 if (curFp >= curFace.size())
1462 curFp = curFace.size()-1;
1467 if (curNb >= nbFace.size())
1473 curNb = nbFace.size()-1;
1476 if (curFace[curFp] != nbFace[curNb])
1480 setPtr->insert(facei);
1481 setPtr->insert(nbFacei);
1512 label nBaffleFaces = 0;
1513 label nErrorDuplicate = 0;
1514 label nErrorOrder = 0;
1517 for (label facei = 0; facei < nFaces(); facei++)
1519 const face& curFace = faces()[facei];
1523 nCommonPoints.
clear();
1527 label pointi = curFace[fp];
1533 label nbFacei = nbs[nbI];
1535 if (facei < nbFacei)
1538 ++(nCommonPoints(nbFacei, 0));
1546 if (checkDuplicateFaces(facei, nCommonPoints, nBaffleFaces, setPtr))
1552 if (checkCommonOrder(facei, nCommonPoints, setPtr))
1564 Info<<
" Number of identical duplicate faces (baffle faces): "
1565 << nBaffleFaces <<
endl;
1568 if (nErrorDuplicate > 0 || nErrorOrder > 0)
1571 if (nErrorDuplicate > 0)
1573 Info<<
" <<Number of duplicate (not baffle) faces found: "
1575 <<
". This might indicate a problem." <<
endl;
1578 if (nErrorOrder > 0)
1580 Info<<
" <<Number of faces with non-consecutive shared points: "
1581 << nErrorOrder <<
". This might indicate a problem." <<
endl;
1587 if (
debug || report)
1589 Info<<
" Face-face connectivity OK." <<
endl;
1612 return checkClosedCells
1630 return checkFaceAreas
1646 return checkCellVolumes
1662 return checkFaceOrthogonality
1675 const scalar minPyrVol,
1679 return checkFacePyramids
1697 return checkFaceSkewness
1712 const scalar maxDeg,
1716 return checkFaceAngles
1730 const scalar warnFlatness,
1734 return checkFaceFlatness
1752 return checkConcaveCells
1764 label nFailedChecks = 0;
1766 if (checkPoints(report)) ++nFailedChecks;
1767 if (checkUpperTriangular(report)) ++nFailedChecks;
1768 if (checkCellsZipUp(report)) ++nFailedChecks;
1769 if (checkFaceVertices(report)) ++nFailedChecks;
1770 if (checkFaceFaces(report)) ++nFailedChecks;
1774 if (
debug || report)
1776 Info<<
" Failed " << nFailedChecks
1777 <<
" mesh topology checks." <<
endl;
1783 if (
debug || report)
1794 label nFailedChecks = 0;
1796 if (checkClosedBoundary(report)) ++nFailedChecks;
1797 if (checkClosedCells(report)) ++nFailedChecks;
1798 if (checkFaceAreas(report)) ++nFailedChecks;
1799 if (checkCellVolumes(report)) ++nFailedChecks;
1800 if (checkFaceOrthogonality(report)) ++nFailedChecks;
1801 if (checkFacePyramids(report)) ++nFailedChecks;
1802 if (checkFaceSkewness(report)) ++nFailedChecks;
1806 if (debug || report)
1808 Info<<
" Failed " << nFailedChecks
1809 <<
" mesh geometry checks." <<
endl;
1815 if (
debug || report)
1832 if (debug || report)
1834 Info<<
" Failed " << nFailedChecks
1835 <<
" mesh checks." <<
endl;
1841 if (
debug || report)
1879 scalar prev = skewThreshold_;
1880 skewThreshold_ = val;
Various functions to operate on Lists.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
void clear()
Remove all entries from table.
void clear()
Clear the list, i.e. set size to zero.
A HashTable to objects of type <T> with a label key.
label size() const noexcept
Number of entries.
A list that is sorted upon construction or when explicitly requested with the sort() method.
const labelList & indices() const noexcept
Return the list of sorted indices. Updated every sort.
void sort()
Forward (stable) sort the list (if changed after construction).
bool empty() const noexcept
True if List is empty (ie, size() is zero).
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
void size(const label n)
Older name for setAddressableSize.
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
static void reduceOr(bool &value, const int communicator=worldComm)
Logical (or) reduction (MPI_AllReduce).
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
A cell is defined as a list of faces with extra functionality.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Class to handle errors and exceptions in a simple, consistent stream-based manner.
A face is a list of labels corresponding to mesh vertices.
label find(const Foam::edge &e) const
Find the edge within the face.
Smooth ATC in cells having a point to a set of patches supplied by type.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
bool checkFaceOrthogonality(const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const
Check for non-orthogonality.
static scalar planarCosAngle_
Threshold where faces are considered coplanar.
virtual const labelList & faceOwner() const =0
Face face-owner addressing.
const labelListList & cellEdges() const
virtual const labelList & faceNeighbour() const =0
Face face-neighbour addressing.
const vectorField & faceCentres() const
virtual bool checkCellsZipUp(const bool report=false, labelHashSet *setPtr=nullptr) const
Check cell zip-up.
static scalar skewThreshold_
Skewness warning threshold.
static scalar setNonOrthThreshold(const scalar)
Set the non-orthogonality warning threshold in degrees.
const labelListList & pointCells() const
const scalarField & cellVolumes() const
virtual const faceList & faces() const =0
Return faces.
label nInternalFaces() const noexcept
Number of internal faces.
virtual bool checkMesh(const bool report=false) const
Check mesh for correctness. Returns false for no error.
virtual bool checkTopology(const bool report=false) const
Check mesh topology for correctness.
bool checkCommonOrder(const label, const Map< label > &, labelHashSet *) const
Check that shared points are in consecutive order.
virtual bool checkUpperTriangular(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face ordering.
static scalar closedThreshold_
Static data to control mesh checking.
const vectorField & cellCentres() const
bool checkFacePyramids(const pointField &points, const vectorField &ctrs, const bool report, const bool detailedReport, const scalar minPyrVol, labelHashSet *setPtr) const
Check face pyramid volume.
bool checkFaceAngles(const pointField &points, const vectorField &faceAreas, const bool report, const scalar maxDeg, labelHashSet *setPtr) const
Check face angles.
static scalar nonOrthThreshold_
Non-orthogonality warning threshold in deg.
static scalar setSkewThreshold(const scalar)
Set the skewness warning threshold as percentage.
static scalar setAspectThreshold(const scalar)
Set the aspect ratio warning threshold.
bool checkDuplicateFaces(const label, const Map< label > &, label &nBaffleFaces, labelHashSet *) const
Check if all points on face are shared with another face.
virtual bool checkFaceVertices(const bool report=false, labelHashSet *setPtr=nullptr) const
Check uniqueness of face vertices.
bool checkCellVolumes(const scalarField &vols, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative cell volumes.
bool checkFaceAreas(const vectorField &faceAreas, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative face areas.
label nFaces() const noexcept
Number of mesh faces.
virtual bool checkPoints(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for unused points.
bool checkFaceSkewness(const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const
Check face skewness.
bool checkFaceFlatness(const pointField &points, const vectorField &faceCentres, const vectorField &faceAreas, const bool report, const scalar warnFlatness, labelHashSet *setPtr) const
Check face warpage.
const labelListList & pointFaces() const
static scalar aspectThreshold_
Aspect ratio warning threshold.
const vectorField & faceAreas() const
bool checkConcaveCells(const vectorField &fAreas, const pointField &fCentres, const bool report, labelHashSet *setPtr) const
Check for concave cells by the planes of faces.
bool checkClosedBoundary(const vectorField &areas, const bool report, const bitSet &internalOrCoupledFaces) const
Check boundary for closedness.
virtual bool checkFaceFaces(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face-face connectivity.
virtual const pointField & points() const =0
Return mesh points.
bool checkClosedCells(const vectorField &faceAreas, const scalarField &cellVolumes, const bool report, labelHashSet *setPtr, labelHashSet *aspectSetPtr, const Vector< label > &meshD) const
Check cells for closedness.
virtual bool checkGeometry(const bool report=false) const
Check mesh geometry (& implicitly topology) for correctness.
static scalar setClosedThreshold(const scalar)
Set the closedness ratio warning threshold.
A class for managing temporary objects.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define DebugInFunction
Report an information message using Foam::Info.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
List< edge > edgeList
List of edge.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar asin(const dimensionedScalar &ds)
Type gSum(const FieldField< Field, Type > &f)
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
label checkGeometry(const polyMesh &mesh, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, autoPtr< coordSetWriter > &setWriter)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
label checkTopology(const polyMesh &mesh, const bool allTopology, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, autoPtr< coordSetWriter > &setWriter, const bool writeBadEdges=false)
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
List< cell > cellList
List of cell.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Vector< solveScalar > solveVector
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Unit conversion functions.