FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

polyMeshGeometry.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "polyMeshGeometry.H"
00027 #include <OpenFOAM/pyramidPointFaceRef.H>
00028 #include <OpenFOAM/syncTools.H>
00029 
00030 namespace Foam
00031 {
00032 
00033 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00034 
00035 defineTypeNameAndDebug(polyMeshGeometry, 0);
00036 
00037 }
00038 
00039 
00040 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00041 
00042 void Foam::polyMeshGeometry::updateFaceCentresAndAreas
00043 (
00044     const pointField& p,
00045     const labelList& changedFaces
00046 )
00047 {
00048     const faceList& fs = mesh_.faces();
00049 
00050     forAll(changedFaces, i)
00051     {
00052         label facei = changedFaces[i];
00053 
00054         const labelList& f = fs[facei];
00055         label nPoints = f.size();
00056 
00057         // If the face is a triangle, do a direct calculation for efficiency
00058         // and to avoid round-off error-related problems
00059         if (nPoints == 3)
00060         {
00061             faceCentres_[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
00062             faceAreas_[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
00063         }
00064         else
00065         {
00066             vector sumN = vector::zero;
00067             scalar sumA = 0.0;
00068             vector sumAc = vector::zero;
00069 
00070             point fCentre = p[f[0]];
00071             for (label pi = 1; pi < nPoints; pi++)
00072             {
00073                 fCentre += p[f[pi]];
00074             }
00075 
00076             fCentre /= nPoints;
00077 
00078             for (label pi = 0; pi < nPoints; pi++)
00079             {
00080                 const point& nextPoint = p[f[(pi + 1) % nPoints]];
00081 
00082                 vector c = p[f[pi]] + nextPoint + fCentre;
00083                 vector n = (nextPoint - p[f[pi]])^(fCentre - p[f[pi]]);
00084                 scalar a = mag(n);
00085 
00086                 sumN += n;
00087                 sumA += a;
00088                 sumAc += a*c;
00089             }
00090 
00091             faceCentres_[facei] = (1.0/3.0)*sumAc/(sumA + VSMALL);
00092             faceAreas_[facei] = 0.5*sumN;
00093         }
00094     }
00095 }
00096 
00097 
00098 void Foam::polyMeshGeometry::updateCellCentresAndVols
00099 (
00100     const labelList& changedCells,
00101     const labelList& changedFaces
00102 )
00103 {
00104     // Clear the fields for accumulation
00105     UIndirectList<vector>(cellCentres_, changedCells) = vector::zero;
00106     UIndirectList<scalar>(cellVolumes_, changedCells) = 0.0;
00107 
00108     const labelList& own = mesh_.faceOwner();
00109     const labelList& nei = mesh_.faceNeighbour();
00110 
00111     // first estimate the approximate cell centre as the average of face centres
00112 
00113     vectorField cEst(mesh_.nCells());
00114     UIndirectList<vector>(cEst, changedCells) = vector::zero;
00115     scalarField nCellFaces(mesh_.nCells());
00116     UIndirectList<scalar>(nCellFaces, changedCells) = 0.0;
00117 
00118     forAll(changedFaces, i)
00119     {
00120         label faceI = changedFaces[i];
00121         cEst[own[faceI]] += faceCentres_[faceI];
00122         nCellFaces[own[faceI]] += 1;
00123 
00124         if (mesh_.isInternalFace(faceI))
00125         {
00126             cEst[nei[faceI]] += faceCentres_[faceI];
00127             nCellFaces[nei[faceI]] += 1;
00128         }
00129     }
00130 
00131     forAll(changedCells, i)
00132     {
00133         label cellI = changedCells[i];
00134         cEst[cellI] /= nCellFaces[cellI];
00135     }
00136 
00137     forAll(changedFaces, i)
00138     {
00139         label faceI = changedFaces[i];
00140 
00141         // Calculate 3*face-pyramid volume
00142         scalar pyr3Vol = max
00143         (
00144             faceAreas_[faceI] & (faceCentres_[faceI] - cEst[own[faceI]]),
00145             VSMALL
00146         );
00147 
00148         // Calculate face-pyramid centre
00149         vector pc = (3.0/4.0)*faceCentres_[faceI] + (1.0/4.0)*cEst[own[faceI]];
00150 
00151         // Accumulate volume-weighted face-pyramid centre
00152         cellCentres_[own[faceI]] += pyr3Vol*pc;
00153 
00154         // Accumulate face-pyramid volume
00155         cellVolumes_[own[faceI]] += pyr3Vol;
00156 
00157         if (mesh_.isInternalFace(faceI))
00158         {
00159             // Calculate 3*face-pyramid volume
00160             scalar pyr3Vol = max
00161             (
00162                 faceAreas_[faceI] & (cEst[nei[faceI]] - faceCentres_[faceI]),
00163                 VSMALL
00164             );
00165 
00166             // Calculate face-pyramid centre
00167             vector pc =
00168                 (3.0/4.0)*faceCentres_[faceI]
00169               + (1.0/4.0)*cEst[nei[faceI]];
00170 
00171             // Accumulate volume-weighted face-pyramid centre
00172             cellCentres_[nei[faceI]] += pyr3Vol*pc;
00173 
00174             // Accumulate face-pyramid volume
00175             cellVolumes_[nei[faceI]] += pyr3Vol;
00176         }
00177     }
00178 
00179     forAll(changedCells, i)
00180     {
00181         label cellI = changedCells[i];
00182 
00183         cellCentres_[cellI] /= cellVolumes_[cellI] + VSMALL;
00184         cellVolumes_[cellI] *= (1.0/3.0);
00185     }
00186 }
00187 
00188 
00189 Foam::labelList Foam::polyMeshGeometry::affectedCells
00190 (
00191     const polyMesh& mesh,
00192     const labelList& changedFaces
00193 )
00194 {
00195     const labelList& own = mesh.faceOwner();
00196     const labelList& nei = mesh.faceNeighbour();
00197 
00198     labelHashSet affectedCells(2*changedFaces.size());
00199 
00200     forAll(changedFaces, i)
00201     {
00202         label faceI = changedFaces[i];
00203 
00204         affectedCells.insert(own[faceI]);
00205 
00206         if (mesh.isInternalFace(faceI))
00207         {
00208             affectedCells.insert(nei[faceI]);
00209         }
00210     }
00211     return affectedCells.toc();
00212 }
00213 
00214 
00215 Foam::scalar Foam::polyMeshGeometry::checkNonOrtho
00216 (
00217     const polyMesh& mesh,
00218     const bool report,
00219     const scalar severeNonorthogonalityThreshold,
00220     const label faceI,
00221     const vector& s,    // face area vector
00222     const vector& d,    // cc-cc vector
00223 
00224     label& severeNonOrth,
00225     label& errorNonOrth,
00226     labelHashSet* setPtr
00227 )
00228 {
00229     scalar dDotS = (d & s)/(mag(d)*mag(s) + VSMALL);
00230 
00231     if (dDotS < severeNonorthogonalityThreshold)
00232     {
00233         label nei = -1;
00234 
00235         if (mesh.isInternalFace(faceI))
00236         {
00237             nei = mesh.faceNeighbour()[faceI];
00238         }
00239 
00240         if (dDotS > SMALL)
00241         {
00242             if (report)
00243             {
00244                 // Severe non-orthogonality but mesh still OK
00245                 Pout<< "Severe non-orthogonality for face " << faceI
00246                     << " between cells " << mesh.faceOwner()[faceI]
00247                     << " and " << nei
00248                     << ": Angle = "
00249                     << ::acos(dDotS)/mathematicalConstant::pi*180.0
00250                     << " deg." << endl;
00251             }
00252 
00253             severeNonOrth++;
00254         }
00255         else
00256         {
00257             // Non-orthogonality greater than 90 deg
00258             if (report)
00259             {
00260                 WarningIn
00261                 (
00262                     "polyMeshGeometry::checkFaceDotProduct"
00263                     "(const bool, const scalar, const labelList&"
00264                     ", labelHashSet*)"
00265                 )   << "Severe non-orthogonality detected for face "
00266                     << faceI
00267                     << " between cells " << mesh.faceOwner()[faceI]
00268                     << " and " << nei
00269                     << ": Angle = "
00270                     << ::acos(dDotS)/mathematicalConstant::pi*180.0
00271                     << " deg." << endl;
00272             }
00273 
00274             errorNonOrth++;
00275         }
00276 
00277         if (setPtr)
00278         {
00279             setPtr->insert(faceI);
00280         }
00281     }
00282     return dDotS;
00283 }
00284 
00285 
00286 Foam::scalar Foam::polyMeshGeometry::calcSkewness
00287 (
00288     const point& ownCc,
00289     const point& neiCc,
00290     const point& fc
00291 )
00292 {
00293     scalar dOwn = mag(fc - ownCc);
00294     scalar dNei = mag(fc - neiCc);
00295 
00296     point faceIntersection =
00297         ownCc*dNei/(dOwn+dNei+VSMALL)
00298       + neiCc*dOwn/(dOwn+dNei+VSMALL);
00299 
00300     return
00301         mag(fc - faceIntersection)
00302       / (
00303             mag(neiCc-ownCc)
00304           + VSMALL
00305         );
00306 }
00307 
00308 
00309 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00310 
00311 // Construct from components
00312 Foam::polyMeshGeometry::polyMeshGeometry(const polyMesh& mesh)
00313 :
00314     mesh_(mesh)
00315 {
00316     correct();
00317 }
00318 
00319 
00320 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00321 
00322 
00323 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00324 
00325 //- Take over properties from mesh
00326 void Foam::polyMeshGeometry::correct()
00327 {
00328     faceAreas_ = mesh_.faceAreas();
00329     faceCentres_ = mesh_.faceCentres();
00330     cellCentres_ = mesh_.cellCentres();
00331     cellVolumes_ = mesh_.cellVolumes();
00332 }
00333 
00334 
00335 //- Recalculate on selected faces
00336 void Foam::polyMeshGeometry::correct
00337 (
00338     const pointField& p,
00339     const labelList& changedFaces
00340 )
00341 {
00342     // Update face quantities
00343     updateFaceCentresAndAreas(p, changedFaces);
00344     // Update cell quantities from face quantities
00345     updateCellCentresAndVols(affectedCells(mesh_, changedFaces), changedFaces);
00346 }
00347 
00348 
00349 bool Foam::polyMeshGeometry::checkFaceDotProduct
00350 (
00351     const bool report,
00352     const scalar orthWarn,
00353     const polyMesh& mesh,
00354     const vectorField& cellCentres,
00355     const vectorField& faceAreas,
00356     const labelList& checkFaces,
00357     const List<labelPair>& baffles,
00358     labelHashSet* setPtr
00359 )
00360 {
00361     // for all internal and coupled faces check theat the d dot S product
00362     // is positive
00363 
00364     const labelList& own = mesh.faceOwner();
00365     const labelList& nei = mesh.faceNeighbour();
00366     const polyBoundaryMesh& patches = mesh.boundaryMesh();
00367 
00368     // Severe nonorthogonality threshold
00369     const scalar severeNonorthogonalityThreshold =
00370         ::cos(orthWarn/180.0*mathematicalConstant::pi);
00371 
00372 
00373     // Calculate coupled cell centre
00374     vectorField neiCc(mesh.nFaces()-mesh.nInternalFaces());
00375 
00376     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
00377     {
00378         neiCc[faceI-mesh.nInternalFaces()] = cellCentres[own[faceI]];
00379     }
00380     syncTools::swapBoundaryFaceList(mesh, neiCc, true);
00381 
00382 
00383     scalar minDDotS = GREAT;
00384 
00385     scalar sumDDotS = 0;
00386     label nDDotS = 0;
00387 
00388     label severeNonOrth = 0;
00389 
00390     label errorNonOrth = 0;
00391 
00392     forAll(checkFaces, i)
00393     {
00394         label faceI = checkFaces[i];
00395 
00396         const point& ownCc = cellCentres[own[faceI]];
00397 
00398         if (mesh.isInternalFace(faceI))
00399         {
00400             scalar dDotS = checkNonOrtho
00401             (
00402                 mesh,
00403                 report,
00404                 severeNonorthogonalityThreshold,
00405                 faceI,
00406                 faceAreas[faceI],
00407                 cellCentres[nei[faceI]] - ownCc,
00408 
00409                 severeNonOrth,
00410                 errorNonOrth,
00411                 setPtr
00412             );
00413 
00414             if (dDotS < minDDotS)
00415             {
00416                 minDDotS = dDotS;
00417             }
00418 
00419             sumDDotS += dDotS;
00420             nDDotS++;
00421         }
00422         else
00423         {
00424             label patchI = patches.whichPatch(faceI);
00425 
00426             if (patches[patchI].coupled())
00427             {
00428                 scalar dDotS = checkNonOrtho
00429                 (
00430                     mesh,
00431                     report,
00432                     severeNonorthogonalityThreshold,
00433                     faceI,
00434                     faceAreas[faceI],
00435                     neiCc[faceI-mesh.nInternalFaces()] - ownCc,
00436 
00437                     severeNonOrth,
00438                     errorNonOrth,
00439                     setPtr
00440                 );
00441 
00442                 if (dDotS < minDDotS)
00443                 {
00444                     minDDotS = dDotS;
00445                 }
00446 
00447                 sumDDotS += dDotS;
00448                 nDDotS++;
00449             }
00450         }
00451     }
00452 
00453     forAll(baffles, i)
00454     {
00455         label face0 = baffles[i].first();
00456         label face1 = baffles[i].second();
00457 
00458         const point& ownCc = cellCentres[own[face0]];
00459 
00460         scalar dDotS = checkNonOrtho
00461         (
00462             mesh,
00463             report,
00464             severeNonorthogonalityThreshold,
00465             face0,
00466             faceAreas[face0],
00467             cellCentres[own[face1]] - ownCc,
00468 
00469             severeNonOrth,
00470             errorNonOrth,
00471             setPtr
00472          );
00473 
00474         if (dDotS < minDDotS)
00475         {
00476             minDDotS = dDotS;
00477         }
00478 
00479         sumDDotS += dDotS;
00480         nDDotS++;
00481     }
00482 
00483     reduce(minDDotS, minOp<scalar>());
00484     reduce(sumDDotS, sumOp<scalar>());
00485     reduce(nDDotS, sumOp<label>());
00486     reduce(severeNonOrth, sumOp<label>());
00487     reduce(errorNonOrth, sumOp<label>());
00488 
00489     // Only report if there are some internal faces
00490     if (nDDotS > 0)
00491     {
00492         if (report && minDDotS < severeNonorthogonalityThreshold)
00493         {
00494             Info<< "Number of non-orthogonality errors: " << errorNonOrth
00495                 << ". Number of severely non-orthogonal faces: "
00496                 << severeNonOrth  << "." << endl;
00497         }
00498     }
00499 
00500     if (report)
00501     {
00502         if (nDDotS > 0)
00503         {
00504             Info<< "Mesh non-orthogonality Max: "
00505                 << ::acos(minDDotS)/mathematicalConstant::pi*180.0
00506                 << " average: " <<
00507                    ::acos(sumDDotS/nDDotS)/mathematicalConstant::pi*180.0
00508                 << endl;
00509         }
00510     }
00511 
00512     if (errorNonOrth > 0)
00513     {
00514         if (report)
00515         {
00516             SeriousErrorIn
00517             (
00518                 "polyMeshGeometry::checkFaceDotProduct"
00519                 "(const bool, const scalar, const labelList&, labelHashSet*)"
00520             )   << "Error in non-orthogonality detected" << endl;
00521         }
00522 
00523         return true;
00524     }
00525     else
00526     {
00527         if (report)
00528         {
00529             Info<< "Non-orthogonality check OK.\n" << endl;
00530         }
00531 
00532         return false;
00533     }
00534 }
00535 
00536 
00537 bool Foam::polyMeshGeometry::checkFacePyramids
00538 (
00539     const bool report,
00540     const scalar minPyrVol,
00541     const polyMesh& mesh,
00542     const vectorField& cellCentres,
00543     const pointField& p,
00544     const labelList& checkFaces,
00545     const List<labelPair>& baffles,
00546     labelHashSet* setPtr
00547 )
00548 {
00549     // check whether face area vector points to the cell with higher label
00550     const labelList& own = mesh.faceOwner();
00551     const labelList& nei = mesh.faceNeighbour();
00552 
00553     const faceList& f = mesh.faces();
00554 
00555     label nErrorPyrs = 0;
00556 
00557     forAll(checkFaces, i)
00558     {
00559         label faceI = checkFaces[i];
00560 
00561         // Create the owner pyramid - it will have negative volume
00562         scalar pyrVol = pyramidPointFaceRef
00563         (
00564             f[faceI],
00565             cellCentres[own[faceI]]
00566         ).mag(p);
00567 
00568         if (pyrVol > -minPyrVol)
00569         {
00570             if (report)
00571             {
00572                 Pout<< "bool polyMeshGeometry::checkFacePyramids("
00573                     << "const bool, const scalar, const pointField&"
00574                     << ", const labelList&, labelHashSet*): "
00575                     << "face " << faceI << " points the wrong way. " << endl
00576                     << "Pyramid volume: " << -pyrVol
00577                     << " Face " << f[faceI] << " area: " << f[faceI].mag(p)
00578                     << " Owner cell: " << own[faceI] << endl
00579                     << "Owner cell vertex labels: "
00580                     << mesh.cells()[own[faceI]].labels(f)
00581                     << endl;
00582             }
00583 
00584 
00585             if (setPtr)
00586             {
00587                 setPtr->insert(faceI);
00588             }
00589 
00590             nErrorPyrs++;
00591         }
00592 
00593         if (mesh.isInternalFace(faceI))
00594         {
00595             // Create the neighbour pyramid - it will have positive volume
00596             scalar pyrVol =
00597                 pyramidPointFaceRef(f[faceI], cellCentres[nei[faceI]]).mag(p);
00598 
00599             if (pyrVol < minPyrVol)
00600             {
00601                 if (report)
00602                 {
00603                     Pout<< "bool polyMeshGeometry::checkFacePyramids("
00604                         << "const bool, const scalar, const pointField&"
00605                         << ", const labelList&, labelHashSet*): "
00606                         << "face " << faceI << " points the wrong way. " << endl
00607                         << "Pyramid volume: " << -pyrVol
00608                         << " Face " << f[faceI] << " area: " << f[faceI].mag(p)
00609                         << " Neighbour cell: " << nei[faceI] << endl
00610                         << "Neighbour cell vertex labels: "
00611                         << mesh.cells()[nei[faceI]].labels(f)
00612                         << endl;
00613                 }
00614 
00615                 if (setPtr)
00616                 {
00617                     setPtr->insert(faceI);
00618                 }
00619 
00620                 nErrorPyrs++;
00621             }
00622         }
00623     }
00624 
00625     forAll(baffles, i)
00626     {
00627         label face0 = baffles[i].first();
00628         label face1 = baffles[i].second();
00629 
00630         const point& ownCc = cellCentres[own[face0]];
00631 
00632        // Create the owner pyramid - it will have negative volume
00633         scalar pyrVolOwn = pyramidPointFaceRef
00634         (
00635             f[face0],
00636             ownCc
00637         ).mag(p);
00638 
00639         if (pyrVolOwn > -minPyrVol)
00640         {
00641             if (report)
00642             {
00643                 Pout<< "bool polyMeshGeometry::checkFacePyramids("
00644                     << "const bool, const scalar, const pointField&"
00645                     << ", const labelList&, labelHashSet*): "
00646                     << "face " << face0 << " points the wrong way. " << endl
00647                     << "Pyramid volume: " << -pyrVolOwn
00648                     << " Face " << f[face0] << " area: " << f[face0].mag(p)
00649                     << " Owner cell: " << own[face0] << endl
00650                     << "Owner cell vertex labels: "
00651                     << mesh.cells()[own[face0]].labels(f)
00652                     << endl;
00653             }
00654 
00655 
00656             if (setPtr)
00657             {
00658                 setPtr->insert(face0);
00659             }
00660 
00661             nErrorPyrs++;
00662         }
00663 
00664         // Create the neighbour pyramid - it will have positive volume
00665         scalar pyrVolNbr =
00666             pyramidPointFaceRef(f[face0], cellCentres[own[face1]]).mag(p);
00667 
00668         if (pyrVolNbr < minPyrVol)
00669         {
00670             if (report)
00671             {
00672                 Pout<< "bool polyMeshGeometry::checkFacePyramids("
00673                     << "const bool, const scalar, const pointField&"
00674                     << ", const labelList&, labelHashSet*): "
00675                     << "face " << face0 << " points the wrong way. " << endl
00676                     << "Pyramid volume: " << -pyrVolNbr
00677                     << " Face " << f[face0] << " area: " << f[face0].mag(p)
00678                     << " Neighbour cell: " << own[face1] << endl
00679                     << "Neighbour cell vertex labels: "
00680                     << mesh.cells()[own[face1]].labels(f)
00681                     << endl;
00682             }
00683 
00684             if (setPtr)
00685             {
00686                 setPtr->insert(face0);
00687             }
00688 
00689             nErrorPyrs++;
00690         }
00691     }
00692 
00693     reduce(nErrorPyrs, sumOp<label>());
00694 
00695     if (nErrorPyrs > 0)
00696     {
00697         if (report)
00698         {
00699             SeriousErrorIn
00700             (
00701                 "polyMeshGeometry::checkFacePyramids("
00702                 "const bool, const scalar, const pointField&"
00703                 ", const labelList&, labelHashSet*)"
00704             )   << "Error in face pyramids: faces pointing the wrong way!"
00705                 << endl;
00706         }
00707 
00708         return true;
00709     }
00710     else
00711     {
00712         if (report)
00713         {
00714             Info<< "Face pyramids OK.\n" << endl;
00715         }
00716 
00717         return false;
00718     }
00719 }
00720 
00721 
00722 bool Foam::polyMeshGeometry::checkFaceSkewness
00723 (
00724     const bool report,
00725     const scalar internalSkew,
00726     const scalar boundarySkew,
00727     const polyMesh& mesh,
00728     const vectorField& cellCentres,
00729     const vectorField& faceCentres,
00730     const vectorField& faceAreas,
00731     const labelList& checkFaces,
00732     const List<labelPair>& baffles,
00733     labelHashSet* setPtr
00734 )
00735 {
00736     // Warn if the skew correction vector is more than skew times
00737     // larger than the face area vector
00738 
00739     const labelList& own = mesh.faceOwner();
00740     const labelList& nei = mesh.faceNeighbour();
00741     const polyBoundaryMesh& patches = mesh.boundaryMesh();
00742 
00743     // Calculate coupled cell centre
00744     vectorField neiCc(mesh.nFaces()-mesh.nInternalFaces());
00745 
00746     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
00747     {
00748         neiCc[faceI-mesh.nInternalFaces()] = cellCentres[own[faceI]];
00749     }
00750     syncTools::swapBoundaryFaceList(mesh, neiCc, true);
00751 
00752 
00753     scalar maxSkew = 0;
00754 
00755     label nWarnSkew = 0;
00756 
00757     forAll(checkFaces, i)
00758     {
00759         label faceI = checkFaces[i];
00760 
00761         if (mesh.isInternalFace(faceI))
00762         {
00763             scalar skewness = calcSkewness
00764             (
00765                 cellCentres[own[faceI]],
00766                 cellCentres[nei[faceI]],
00767                 faceCentres[faceI]
00768             );
00769 
00770             // Check if the skewness vector is greater than the PN vector.
00771             // This does not cause trouble but is a good indication of a poor
00772             // mesh.
00773             if (skewness > internalSkew)
00774             {
00775                 if (report)
00776                 {
00777                     Pout<< "Severe skewness for face " << faceI
00778                         << " skewness = " << skewness << endl;
00779                 }
00780 
00781                 if (setPtr)
00782                 {
00783                     setPtr->insert(faceI);
00784                 }
00785 
00786                 nWarnSkew++;
00787             }
00788 
00789             maxSkew = max(maxSkew, skewness);
00790         }
00791         else if (patches[patches.whichPatch(faceI)].coupled())
00792         {
00793             scalar skewness = calcSkewness
00794             (
00795                 cellCentres[own[faceI]],
00796                 neiCc[faceI-mesh.nInternalFaces()],
00797                 faceCentres[faceI]
00798             );
00799 
00800             // Check if the skewness vector is greater than the PN vector.
00801             // This does not cause trouble but is a good indication of a poor
00802             // mesh.
00803             if (skewness > internalSkew)
00804             {
00805                 if (report)
00806                 {
00807                     Pout<< "Severe skewness for coupled face " << faceI
00808                         << " skewness = " << skewness << endl;
00809                 }
00810 
00811                 if (setPtr)
00812                 {
00813                     setPtr->insert(faceI);
00814                 }
00815 
00816                 nWarnSkew++;
00817             }
00818 
00819             maxSkew = max(maxSkew, skewness);
00820         }
00821         else
00822         {
00823             // Boundary faces: consider them to have only skewness error.
00824             // (i.e. treat as if mirror cell on other side)
00825 
00826             vector faceNormal = faceAreas[faceI];
00827             faceNormal /= mag(faceNormal) + ROOTVSMALL;
00828 
00829             vector dOwn = faceCentres[faceI] - cellCentres[own[faceI]];
00830 
00831             vector dWall = faceNormal*(faceNormal & dOwn);
00832 
00833             point faceIntersection = cellCentres[own[faceI]] + dWall;
00834 
00835             scalar skewness =
00836                 mag(faceCentres[faceI] - faceIntersection)
00837                 /(2*mag(dWall) + ROOTVSMALL);
00838 
00839             // Check if the skewness vector is greater than the PN vector.
00840             // This does not cause trouble but is a good indication of a poor
00841             // mesh.
00842             if (skewness > boundarySkew)
00843             {
00844                 if (report)
00845                 {
00846                     Pout<< "Severe skewness for boundary face " << faceI
00847                         << " skewness = " << skewness << endl;
00848                 }
00849 
00850                 if (setPtr)
00851                 {
00852                     setPtr->insert(faceI);
00853                 }
00854 
00855                 nWarnSkew++;
00856             }
00857 
00858             maxSkew = max(maxSkew, skewness);
00859         }
00860     }
00861 
00862     forAll(baffles, i)
00863     {
00864         label face0 = baffles[i].first();
00865         label face1 = baffles[i].second();
00866 
00867         const point& ownCc = cellCentres[own[face0]];
00868 
00869         scalar skewness = calcSkewness
00870         (
00871             ownCc,
00872             cellCentres[own[face1]],
00873             faceCentres[face0]
00874         );
00875 
00876         // Check if the skewness vector is greater than the PN vector.
00877         // This does not cause trouble but is a good indication of a poor
00878         // mesh.
00879         if (skewness > internalSkew)
00880         {
00881             if (report)
00882             {
00883                 Pout<< "Severe skewness for face " << face0
00884                     << " skewness = " << skewness << endl;
00885             }
00886 
00887             if (setPtr)
00888             {
00889                 setPtr->insert(face0);
00890             }
00891 
00892             nWarnSkew++;
00893         }
00894 
00895         maxSkew = max(maxSkew, skewness);
00896     }
00897 
00898 
00899     reduce(maxSkew, maxOp<scalar>());
00900     reduce(nWarnSkew, sumOp<label>());
00901 
00902     if (nWarnSkew > 0)
00903     {
00904         if (report)
00905         {
00906             WarningIn
00907             (
00908                 "polyMeshGeometry::checkFaceSkewness"
00909                 "(const bool, const scalar, const labelList&, labelHashSet*)"
00910             )   << "Large face skewness detected.  Max skewness = "
00911                 << 100*maxSkew
00912                 << " percent.\nThis may impair the quality of the result." << nl
00913                 << nWarnSkew << " highly skew faces detected."
00914                 << endl;
00915         }
00916 
00917         return true;
00918     }
00919     else
00920     {
00921         if (report)
00922         {
00923             Info<< "Max skewness = " << 100*maxSkew
00924                 << " percent.  Face skewness OK.\n" << endl;
00925         }
00926 
00927         return false;
00928     }
00929 }
00930 
00931 
00932 bool Foam::polyMeshGeometry::checkFaceWeights
00933 (
00934     const bool report,
00935     const scalar warnWeight,
00936     const polyMesh& mesh,
00937     const vectorField& cellCentres,
00938     const vectorField& faceCentres,
00939     const vectorField& faceAreas,
00940     const labelList& checkFaces,
00941     const List<labelPair>& baffles,
00942     labelHashSet* setPtr
00943 )
00944 {
00945     // Warn if the delta factor (0..1) is too large.
00946 
00947     const labelList& own = mesh.faceOwner();
00948     const labelList& nei = mesh.faceNeighbour();
00949     const polyBoundaryMesh& patches = mesh.boundaryMesh();
00950 
00951     // Calculate coupled cell centre
00952     vectorField neiCc(mesh.nFaces()-mesh.nInternalFaces());
00953 
00954     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
00955     {
00956         neiCc[faceI-mesh.nInternalFaces()] = cellCentres[own[faceI]];
00957     }
00958     syncTools::swapBoundaryFaceList(mesh, neiCc, true);
00959 
00960 
00961     scalar minWeight = GREAT;
00962 
00963     label nWarnWeight = 0;
00964 
00965     forAll(checkFaces, i)
00966     {
00967         label faceI = checkFaces[i];
00968 
00969         const point& fc = faceCentres[faceI];
00970         const vector& fa = faceAreas[faceI];
00971 
00972         scalar dOwn = mag(fa & (fc-cellCentres[own[faceI]]));
00973 
00974         if (mesh.isInternalFace(faceI))
00975         {
00976             scalar dNei = mag(fa & (cellCentres[nei[faceI]]-fc));
00977             scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
00978 
00979             if (weight < warnWeight)
00980             {
00981                 if (report)
00982                 {
00983                     Pout<< "Small weighting factor for face " << faceI
00984                         << " weight = " << weight << endl;
00985                 }
00986 
00987                 if (setPtr)
00988                 {
00989                     setPtr->insert(faceI);
00990                 }
00991 
00992                 nWarnWeight++;
00993             }
00994 
00995             minWeight = min(minWeight, weight);
00996         }
00997         else
00998         {
00999             label patchI = patches.whichPatch(faceI);
01000 
01001             if (patches[patchI].coupled())
01002             {
01003                 scalar dNei = mag(fa & (neiCc[faceI-mesh.nInternalFaces()]-fc));
01004                 scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
01005 
01006                 if (weight < warnWeight)
01007                 {
01008                     if (report)
01009                     {
01010                         Pout<< "Small weighting factor for face " << faceI
01011                             << " weight = " << weight << endl;
01012                     }
01013 
01014                     if (setPtr)
01015                     {
01016                         setPtr->insert(faceI);
01017                     }
01018 
01019                     nWarnWeight++;
01020                 }
01021 
01022                 minWeight = min(minWeight, weight);
01023             }
01024         }
01025     }
01026 
01027     forAll(baffles, i)
01028     {
01029         label face0 = baffles[i].first();
01030         label face1 = baffles[i].second();
01031 
01032         const point& ownCc = cellCentres[own[face0]];
01033         const point& fc = faceCentres[face0];
01034         const vector& fa = faceAreas[face0];
01035 
01036         scalar dOwn = mag(fa & (fc-ownCc));
01037         scalar dNei = mag(fa & (cellCentres[own[face1]]-fc));
01038         scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
01039 
01040         if (weight < warnWeight)
01041         {
01042             if (report)
01043             {
01044                 Pout<< "Small weighting factor for face " << face0
01045                     << " weight = " << weight << endl;
01046             }
01047 
01048             if (setPtr)
01049             {
01050                 setPtr->insert(face0);
01051             }
01052 
01053             nWarnWeight++;
01054         }
01055 
01056         minWeight = min(minWeight, weight);
01057     }
01058 
01059     reduce(minWeight, minOp<scalar>());
01060     reduce(nWarnWeight, sumOp<label>());
01061 
01062     if (minWeight < warnWeight)
01063     {
01064         if (report)
01065         {
01066             WarningIn
01067             (
01068                 "polyMeshGeometry::checkFaceWeights"
01069                 "(const bool, const scalar, const labelList&, labelHashSet*)"
01070             )   << "Small interpolation weight detected.  Min weight = "
01071                 << minWeight << '.' << nl
01072                 << nWarnWeight << " faces with small weights detected."
01073                 << endl;
01074         }
01075 
01076         return true;
01077     }
01078     else
01079     {
01080         if (report)
01081         {
01082             Info<< "Min weight = " << minWeight
01083                 << ".  Weights OK.\n" << endl;
01084         }
01085 
01086         return false;
01087     }
01088 }
01089 
01090 
01091 bool Foam::polyMeshGeometry::checkVolRatio
01092 (
01093     const bool report,
01094     const scalar warnRatio,
01095     const polyMesh& mesh,
01096     const scalarField& cellVolumes,
01097     const labelList& checkFaces,
01098     const List<labelPair>& baffles,
01099     labelHashSet* setPtr
01100 )
01101 {
01102     // Warn if the volume ratio between neighbouring cells is too large
01103 
01104     const labelList& own = mesh.faceOwner();
01105     const labelList& nei = mesh.faceNeighbour();
01106     const polyBoundaryMesh& patches = mesh.boundaryMesh();
01107 
01108     // Calculate coupled cell vol
01109     scalarField neiVols(mesh.nFaces()-mesh.nInternalFaces());
01110 
01111     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
01112     {
01113         neiVols[faceI-mesh.nInternalFaces()] = cellVolumes[own[faceI]];
01114     }
01115     syncTools::swapBoundaryFaceList(mesh, neiVols, true);
01116 
01117 
01118     scalar minRatio = GREAT;
01119 
01120     label nWarnRatio = 0;
01121 
01122     forAll(checkFaces, i)
01123     {
01124         label faceI = checkFaces[i];
01125 
01126         scalar ownVol = mag(cellVolumes[own[faceI]]);
01127 
01128         scalar neiVol = -GREAT;
01129 
01130         if (mesh.isInternalFace(faceI))
01131         {
01132             neiVol = mag(cellVolumes[nei[faceI]]);
01133         }
01134         else
01135         {
01136             label patchI = patches.whichPatch(faceI);
01137 
01138             if (patches[patchI].coupled())
01139             {
01140                 neiVol = mag(neiVols[faceI-mesh.nInternalFaces()]);
01141             }
01142         }
01143 
01144         if (neiVol >= 0)
01145         {
01146             scalar ratio = min(ownVol, neiVol) / (max(ownVol, neiVol) + VSMALL);
01147 
01148             if (ratio < warnRatio)
01149             {
01150                 if (report)
01151                 {
01152                     Pout<< "Small ratio for face " << faceI
01153                         << " ratio = " << ratio << endl;
01154                 }
01155 
01156                 if (setPtr)
01157                 {
01158                     setPtr->insert(faceI);
01159                 }
01160 
01161                 nWarnRatio++;
01162             }
01163 
01164             minRatio = min(minRatio, ratio);
01165         }
01166     }
01167 
01168     forAll(baffles, i)
01169     {
01170         label face0 = baffles[i].first();
01171         label face1 = baffles[i].second();
01172         
01173         scalar ownVol = mag(cellVolumes[own[face0]]);
01174 
01175         scalar neiVol = mag(cellVolumes[own[face1]]);
01176 
01177         if (neiVol >= 0)
01178         {
01179             scalar ratio = min(ownVol, neiVol) / (max(ownVol, neiVol) + VSMALL);
01180 
01181             if (ratio < warnRatio)
01182             {
01183                 if (report)
01184                 {
01185                     Pout<< "Small ratio for face " << face0
01186                         << " ratio = " << ratio << endl;
01187                 }
01188 
01189                 if (setPtr)
01190                 {
01191                     setPtr->insert(face0);
01192                 }
01193 
01194                 nWarnRatio++;
01195             }
01196 
01197             minRatio = min(minRatio, ratio);
01198         }
01199     }
01200 
01201     reduce(minRatio, minOp<scalar>());
01202     reduce(nWarnRatio, sumOp<label>());
01203 
01204     if (minRatio < warnRatio)
01205     {
01206         if (report)
01207         {
01208             WarningIn
01209             (
01210                 "polyMeshGeometry::checkVolRatio"
01211                 "(const bool, const scalar, const labelList&, labelHashSet*)"
01212             )   << "Small volume ratio detected.  Min ratio = "
01213                 << minRatio << '.' << nl
01214                 << nWarnRatio << " faces with small ratios detected."
01215                 << endl;
01216         }
01217 
01218         return true;
01219     }
01220     else
01221     {
01222         if (report)
01223         {
01224             Info<< "Min ratio = " << minRatio
01225                 << ".  Ratios OK.\n" << endl;
01226         }
01227 
01228         return false;
01229     }
01230 }
01231 
01232 
01233 // Check convexity of angles in a face. Allow a slight non-convexity.
01234 // E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity)
01235 // (if truly concave and points not visible from face centre the face-pyramid
01236 //  check in checkMesh will fail)
01237 bool Foam::polyMeshGeometry::checkFaceAngles
01238 (
01239     const bool report,
01240     const scalar maxDeg,
01241     const polyMesh& mesh,
01242     const vectorField& faceAreas,
01243     const pointField& p,
01244     const labelList& checkFaces,
01245     labelHashSet* setPtr
01246 )
01247 {
01248     if (maxDeg < -SMALL || maxDeg > 180+SMALL)
01249     {
01250         FatalErrorIn
01251         (
01252             "polyMeshGeometry::checkFaceAngles"
01253             "(const bool, const scalar, const pointField&, const labelList&"
01254             ", labelHashSet*)"
01255         )   << "maxDeg should be [0..180] but is now " << maxDeg
01256             << abort(FatalError);
01257     }
01258 
01259     const scalar maxSin = Foam::sin(maxDeg/180.0*mathematicalConstant::pi);
01260 
01261     const faceList& fcs = mesh.faces();
01262 
01263     scalar maxEdgeSin = 0.0;
01264 
01265     label nConcave = 0;
01266 
01267     label errorFaceI = -1;
01268 
01269     forAll(checkFaces, i)
01270     {
01271         label faceI = checkFaces[i];
01272 
01273         const face& f = fcs[faceI];
01274 
01275         vector faceNormal = faceAreas[faceI];
01276         faceNormal /= mag(faceNormal) + VSMALL;
01277 
01278         // Get edge from f[0] to f[size-1];
01279         vector ePrev(p[f[0]] - p[f[f.size()-1]]);
01280         scalar magEPrev = mag(ePrev);
01281         ePrev /= magEPrev + VSMALL;
01282 
01283         forAll(f, fp0)
01284         {
01285             // Get vertex after fp
01286             label fp1 = f.fcIndex(fp0);
01287 
01288             // Normalized vector between two consecutive points
01289             vector e10(p[f[fp1]] - p[f[fp0]]);
01290             scalar magE10 = mag(e10);
01291             e10 /= magE10 + VSMALL;
01292 
01293             if (magEPrev > SMALL && magE10 > SMALL)
01294             {
01295                 vector edgeNormal = ePrev ^ e10;
01296                 scalar magEdgeNormal = mag(edgeNormal);
01297 
01298                 if (magEdgeNormal < maxSin)
01299                 {
01300                     // Edges (almost) aligned -> face is ok.
01301                 }
01302                 else
01303                 {
01304                     // Check normal
01305                     edgeNormal /= magEdgeNormal;
01306 
01307                     if ((edgeNormal & faceNormal) < SMALL)
01308                     {
01309                         if (faceI != errorFaceI)
01310                         {
01311                             // Count only one error per face.
01312                             errorFaceI = faceI;
01313                             nConcave++;
01314                         }
01315 
01316                         if (setPtr)
01317                         {
01318                             setPtr->insert(faceI);
01319                         }
01320 
01321                         maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
01322                     }
01323                 }
01324             }
01325 
01326             ePrev = e10;
01327             magEPrev = magE10;
01328         }
01329     }
01330 
01331     reduce(nConcave, sumOp<label>());
01332     reduce(maxEdgeSin, maxOp<scalar>());
01333 
01334     if (report)
01335     {
01336         if (maxEdgeSin > SMALL)
01337         {
01338             scalar maxConcaveDegr =
01339                 Foam::asin(Foam::min(1.0, maxEdgeSin))
01340              * 180.0/mathematicalConstant::pi;
01341 
01342             Info<< "There are " << nConcave
01343                 << " faces with concave angles between consecutive"
01344                 << " edges. Max concave angle = " << maxConcaveDegr
01345                 << " degrees.\n" << endl;
01346         }
01347         else
01348         {
01349             Info<< "All angles in faces are convex or less than "  << maxDeg
01350                 << " degrees concave.\n" << endl;
01351         }
01352     }
01353 
01354     if (nConcave > 0)
01355     {
01356         if (report)
01357         {
01358             WarningIn
01359             (
01360                 "polyMeshGeometry::checkFaceAngles"
01361                 "(const bool, const scalar,  const pointField&"
01362                 ", const labelList&, labelHashSet*)"
01363             )   << nConcave  << " face points with severe concave angle (> "
01364                 << maxDeg << " deg) found.\n"
01365                 << endl;
01366         }
01367 
01368         return true;
01369     }
01370     else
01371     {
01372         return false;
01373     }
01374 }
01375 
01376 
01377 // Check twist of faces. Is calculated as the difference between normals of
01378 // individual triangles and the cell-cell centre edge.
01379 bool Foam::polyMeshGeometry::checkFaceTwist
01380 (
01381     const bool report,
01382     const scalar minTwist,
01383     const polyMesh& mesh,
01384     const vectorField& cellCentres,
01385     const vectorField& faceAreas,
01386     const vectorField& faceCentres,
01387     const pointField& p,
01388     const labelList& checkFaces,
01389     labelHashSet* setPtr
01390 )
01391 {
01392     if (minTwist < -1-SMALL || minTwist > 1+SMALL)
01393     {
01394         FatalErrorIn
01395         (
01396             "polyMeshGeometry::checkFaceTwist"
01397             "(const bool, const scalar, const polyMesh&, const pointField&"
01398             ", const pointField&, const pointField&, const pointField&"
01399             ", const labelList&, labelHashSet*)"
01400         )   << "minTwist should be [-1..1] but is now " << minTwist
01401             << abort(FatalError);
01402     }
01403 
01404 
01405     const faceList& fcs = mesh.faces();
01406 
01407 
01408     label nWarped = 0;
01409 
01410 //    forAll(checkFaces, i)
01411 //    {
01412 //        label faceI = checkFaces[i];
01413 //
01414 //        const face& f = fcs[faceI];
01415 //
01416 //        scalar magArea = mag(faceAreas[faceI]);
01417 //
01418 //        if (f.size() > 3 && magArea > VSMALL)
01419 //        {
01420 //            const vector nf = faceAreas[faceI] / magArea;
01421 //
01422 //            const point& fc = faceCentres[faceI];
01423 //
01424 //            forAll(f, fpI)
01425 //            {
01426 //                vector triArea
01427 //                (
01428 //                    triPointRef
01429 //                    (
01430 //                        p[f[fpI]],
01431 //                        p[f.nextLabel(fpI)],
01432 //                        fc
01433 //                    ).normal()
01434 //                );
01435 //
01436 //                scalar magTri = mag(triArea);
01437 //
01438 //                if (magTri > VSMALL && ((nf & triArea/magTri) < minTwist))
01439 //                {
01440 //                    nWarped++;
01441 //
01442 //                    if (setPtr)
01443 //                    {
01444 //                        setPtr->insert(faceI);
01445 //                    }
01446 //
01447 //                    break;
01448 //                }
01449 //            }
01450 //        }
01451 //    }
01452 
01453 
01454     const labelList& own = mesh.faceOwner();
01455     const labelList& nei = mesh.faceNeighbour();
01456     const polyBoundaryMesh& patches = mesh.boundaryMesh();
01457 
01458     // Calculate coupled cell centre
01459     vectorField neiCc(mesh.nFaces()-mesh.nInternalFaces());
01460 
01461     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
01462     {
01463         neiCc[faceI-mesh.nInternalFaces()] = cellCentres[own[faceI]];
01464     }
01465     syncTools::swapBoundaryFaceList(mesh, neiCc, true);
01466 
01467     forAll(checkFaces, i)
01468     {
01469         label faceI = checkFaces[i];
01470 
01471         const face& f = fcs[faceI];
01472 
01473         if (f.size() > 3)
01474         {
01475             vector nf(vector::zero);
01476 
01477             if (mesh.isInternalFace(faceI))
01478             {
01479                 nf = cellCentres[nei[faceI]] - cellCentres[own[faceI]];
01480                 nf /= mag(nf) + VSMALL;
01481             }
01482             else if (patches[patches.whichPatch(faceI)].coupled())
01483             {
01484                 nf =
01485                     neiCc[faceI-mesh.nInternalFaces()]
01486                   - cellCentres[own[faceI]];
01487                 nf /= mag(nf) + VSMALL;
01488             }
01489             else
01490             {
01491                 nf = faceCentres[faceI] - cellCentres[own[faceI]];
01492                 nf /= mag(nf) + VSMALL;
01493             }
01494 
01495             if (nf != vector::zero)
01496             {
01497                 const point& fc = faceCentres[faceI];
01498 
01499                 forAll(f, fpI)
01500                 {
01501                     vector triArea
01502                     (
01503                         triPointRef
01504                         (
01505                             p[f[fpI]],
01506                             p[f.nextLabel(fpI)],
01507                             fc
01508                         ).normal()
01509                     );
01510 
01511                     scalar magTri = mag(triArea);
01512 
01513                     if (magTri > VSMALL && ((nf & triArea/magTri) < minTwist))
01514                     {
01515                         nWarped++;
01516 
01517                         if (setPtr)
01518                         {
01519                             setPtr->insert(faceI);
01520                         }
01521 
01522                         break;
01523                     }
01524                 }
01525             }
01526         }
01527     }
01528 
01529     reduce(nWarped, sumOp<label>());
01530 
01531     if (report)
01532     {
01533         if (nWarped> 0)
01534         {
01535             Info<< "There are " << nWarped
01536                 << " faces with cosine of the angle"
01537                 << " between triangle normal and face normal less than "
01538                 << minTwist << nl << endl;
01539         }
01540         else
01541         {
01542             Info<< "All faces are flat in that the cosine of the angle"
01543                 << " between triangle normal and face normal less than "
01544                 << minTwist << nl << endl;
01545         }
01546     }
01547 
01548     if (nWarped > 0)
01549     {
01550         if (report)
01551         {
01552             WarningIn
01553             (
01554                 "polyMeshGeometry::checkFaceTwist"
01555                 "(const bool, const scalar, const polyMesh&, const pointField&"
01556                 ", const pointField&, const pointField&, const pointField&"
01557                 ", const labelList&, labelHashSet*)"
01558             )   << nWarped  << " faces with severe warpage "
01559                 << "(cosine of the angle between triangle normal and face normal"
01560                 << " < " << minTwist << ") found.\n"
01561                 << endl;
01562         }
01563 
01564         return true;
01565     }
01566     else
01567     {
01568         return false;
01569     }
01570 }
01571 
01572 
01573 // Like checkFaceTwist but compares normals of consecutive triangles.
01574 bool Foam::polyMeshGeometry::checkTriangleTwist
01575 (
01576     const bool report,
01577     const scalar minTwist,
01578     const polyMesh& mesh,
01579     const vectorField& faceAreas,
01580     const vectorField& faceCentres,
01581     const pointField& p,
01582     const labelList& checkFaces,
01583     labelHashSet* setPtr
01584 )
01585 {
01586     if (minTwist < -1-SMALL || minTwist > 1+SMALL)
01587     {
01588         FatalErrorIn
01589         (
01590             "polyMeshGeometry::checkTriangleTwist"
01591             "(const bool, const scalar, const polyMesh&, const pointField&"
01592             ", const labelList&, labelHashSet*)"
01593         )   << "minTwist should be [-1..1] but is now " << minTwist
01594             << abort(FatalError);
01595     }
01596 
01597     const faceList& fcs = mesh.faces();
01598 
01599     label nWarped = 0;
01600 
01601     forAll(checkFaces, i)
01602     {
01603         label faceI = checkFaces[i];
01604 
01605         const face& f = fcs[faceI];
01606 
01607         if (f.size() > 3)
01608         {
01609             const point& fc = faceCentres[faceI];
01610 
01611             // Find starting triangle (at startFp) with non-zero area
01612             label startFp = -1;
01613             vector prevN;
01614 
01615             forAll(f, fp)
01616             {
01617                 prevN = triPointRef
01618                 (
01619                     p[f[fp]],
01620                     p[f.nextLabel(fp)],
01621                     fc
01622                 ).normal();
01623 
01624                 scalar magTri = mag(prevN);
01625 
01626                 if (magTri > VSMALL)
01627                 {
01628                     startFp = fp;
01629                     prevN /= magTri;
01630                     break;
01631                 }
01632             }
01633 
01634             if (startFp != -1)
01635             {
01636                 label fp = startFp;
01637 
01638                 do
01639                 {
01640                     fp = f.fcIndex(fp);
01641 
01642                     vector triN
01643                     (
01644                         triPointRef
01645                         (
01646                             p[f[fp]],
01647                             p[f.nextLabel(fp)],
01648                             fc
01649                         ).normal()
01650                     );
01651                     scalar magTri = mag(triN);
01652 
01653                     if (magTri > VSMALL)
01654                     {
01655                         triN /= magTri;
01656 
01657                         if ((prevN & triN) < minTwist)
01658                         {
01659                             nWarped++;
01660 
01661                             if (setPtr)
01662                             {
01663                                 setPtr->insert(faceI);
01664                             }
01665 
01666                             break;
01667                         }
01668 
01669                         prevN = triN;
01670                     }
01671                     else if (minTwist > 0)
01672                     {
01673                         nWarped++;
01674 
01675                         if (setPtr)
01676                         {
01677                             setPtr->insert(faceI);
01678                         }
01679 
01680                         break;
01681                     }
01682                 }
01683                 while (fp != startFp);
01684             }
01685         }
01686     }
01687 
01688 
01689     reduce(nWarped, sumOp<label>());
01690 
01691     if (report)
01692     {
01693         if (nWarped> 0)
01694         {
01695             Info<< "There are " << nWarped
01696                 << " faces with cosine of the angle"
01697                 << " between consecutive triangle normals less than "
01698                 << minTwist << nl << endl;
01699         }
01700         else
01701         {
01702             Info<< "All faces are flat in that the cosine of the angle"
01703                 << " between consecutive triangle normals is less than "
01704                 << minTwist << nl << endl;
01705         }
01706     }
01707 
01708     if (nWarped > 0)
01709     {
01710         if (report)
01711         {
01712             WarningIn
01713             (
01714                 "polyMeshGeometry::checkTriangleTwist"
01715                 "(const bool, const scalar, const polyMesh&"
01716                 ", const pointField&, const labelList&, labelHashSet*)"
01717             )   << nWarped  << " faces with severe warpage "
01718                 << "(cosine of the angle between consecutive triangle normals"
01719                 << " < " << minTwist << ") found.\n"
01720                 << endl;
01721         }
01722 
01723         return true;
01724     }
01725     else
01726     {
01727         return false;
01728     }
01729 }
01730 
01731 
01732 bool Foam::polyMeshGeometry::checkFaceArea
01733 (
01734     const bool report,
01735     const scalar minArea,
01736     const polyMesh& mesh,
01737     const vectorField& faceAreas,
01738     const labelList& checkFaces,
01739     labelHashSet* setPtr
01740 )
01741 {
01742     label nZeroArea = 0;
01743 
01744     forAll(checkFaces, i)
01745     {
01746         label faceI = checkFaces[i];
01747 
01748         if (mag(faceAreas[faceI]) < minArea)
01749         {
01750             if (setPtr)
01751             {
01752                 setPtr->insert(faceI);
01753             }
01754             nZeroArea++;
01755         }
01756     }
01757 
01758 
01759     reduce(nZeroArea, sumOp<label>());
01760 
01761     if (report)
01762     {
01763         if (nZeroArea > 0)
01764         {
01765             Info<< "There are " << nZeroArea
01766                 << " faces with area < " << minArea << '.' << nl << endl;
01767         }
01768         else
01769         {
01770             Info<< "All faces have area > " << minArea << '.' << nl << endl;
01771         }
01772     }
01773 
01774     if (nZeroArea > 0)
01775     {
01776         if (report)
01777         {
01778             WarningIn
01779             (
01780                 "polyMeshGeometry::checkFaceArea"
01781                 "(const bool, const scalar, const polyMesh&"
01782                 ", const pointField&, const labelList&, labelHashSet*)"
01783             )   << nZeroArea  << " faces with area < " << minArea
01784                 << " found.\n"
01785                 << endl;
01786         }
01787 
01788         return true;
01789     }
01790     else
01791     {
01792         return false;
01793     }
01794 }
01795 
01796 
01797 bool Foam::polyMeshGeometry::checkCellDeterminant
01798 (
01799     const bool report,
01800     const scalar warnDet,
01801     const polyMesh& mesh,
01802     const vectorField& faceAreas,
01803     const labelList& checkFaces,
01804     const labelList& affectedCells,
01805     labelHashSet* setPtr
01806 )
01807 {
01808     const cellList& cells = mesh.cells();
01809 
01810     scalar minDet = GREAT;
01811     scalar sumDet = 0.0;
01812     label nSumDet = 0;
01813     label nWarnDet = 0;
01814 
01815     forAll(affectedCells, i)
01816     {
01817         const cell& cFaces = cells[affectedCells[i]];
01818 
01819         tensor areaSum(tensor::zero);
01820         scalar magAreaSum = 0;
01821 
01822         forAll(cFaces, cFaceI)
01823         {
01824             label faceI = cFaces[cFaceI];
01825     
01826             scalar magArea = mag(faceAreas[faceI]);
01827 
01828             magAreaSum += magArea;
01829             areaSum += faceAreas[faceI]*(faceAreas[faceI]/(magArea+VSMALL));
01830         }
01831 
01832         scalar scaledDet = det(areaSum/(magAreaSum+VSMALL))/0.037037037037037;
01833 
01834         minDet = min(minDet, scaledDet);
01835         sumDet += scaledDet;
01836         nSumDet++;
01837 
01838         if (scaledDet < warnDet)
01839         {
01840             if (setPtr)
01841             {
01842                 // Insert all faces of the cell.
01843                 forAll(cFaces, cFaceI)
01844                 {
01845                     label faceI = cFaces[cFaceI];
01846                     setPtr->insert(faceI);
01847                 }
01848             }
01849             nWarnDet++;
01850         }
01851     }
01852     
01853     reduce(minDet, minOp<scalar>());
01854     reduce(sumDet, sumOp<scalar>());
01855     reduce(nSumDet, sumOp<label>());
01856     reduce(nWarnDet, sumOp<label>());
01857 
01858     if (report)
01859     {
01860         if (nSumDet > 0)
01861         {
01862             Info<< "Cell determinant (1 = uniform cube) : average = "
01863                 << sumDet / nSumDet << "  min = " << minDet << endl;
01864         }
01865 
01866         if (nWarnDet > 0)
01867         {
01868             Info<< "There are " << nWarnDet
01869                 << " cells with determinant < " << warnDet << '.' << nl
01870                 << endl;
01871         }
01872         else
01873         {
01874             Info<< "All faces have determinant > " << warnDet << '.' << nl
01875                 << endl;
01876         }
01877     }
01878 
01879     if (nWarnDet > 0)
01880     {
01881         if (report)
01882         {
01883             WarningIn
01884             (
01885                 "polyMeshGeometry::checkCellDeterminant"
01886                 "(const bool, const scalar, const polyMesh&"
01887                 ", const pointField&, const labelList&, const labelList&"
01888                 ", labelHashSet*)"
01889             )   << nWarnDet << " cells with determinant < " << warnDet
01890                 << " found.\n"
01891                 << endl;
01892         }
01893 
01894         return true;
01895     }
01896     else
01897     {
01898         return false;
01899     }
01900 }
01901 
01902 
01903 bool Foam::polyMeshGeometry::checkFaceDotProduct
01904 (
01905     const bool report,
01906     const scalar orthWarn,
01907     const labelList& checkFaces,
01908     const List<labelPair>& baffles,
01909     labelHashSet* setPtr
01910 ) const
01911 {
01912     return checkFaceDotProduct
01913     (
01914         report,
01915         orthWarn,
01916         mesh_,
01917         cellCentres_,
01918         faceAreas_,
01919         checkFaces,
01920         baffles,
01921         setPtr
01922     );
01923 }
01924 
01925 
01926 bool Foam::polyMeshGeometry::checkFacePyramids
01927 (
01928     const bool report,
01929     const scalar minPyrVol,
01930     const pointField& p,
01931     const labelList& checkFaces,
01932     const List<labelPair>& baffles,
01933     labelHashSet* setPtr
01934 ) const
01935 {
01936     return checkFacePyramids
01937     (
01938         report,
01939         minPyrVol,
01940         mesh_,
01941         cellCentres_,
01942         p,
01943         checkFaces,
01944         baffles,
01945         setPtr
01946     );
01947 }
01948 
01949 
01950 bool Foam::polyMeshGeometry::checkFaceSkewness
01951 (
01952     const bool report,
01953     const scalar internalSkew,
01954     const scalar boundarySkew,
01955     const labelList& checkFaces,
01956     const List<labelPair>& baffles,
01957     labelHashSet* setPtr
01958 ) const
01959 {
01960     return checkFaceSkewness
01961     (
01962         report,
01963         internalSkew,
01964         boundarySkew,
01965         mesh_,
01966         cellCentres_,
01967         faceCentres_,
01968         faceAreas_,
01969         checkFaces,
01970         baffles,
01971         setPtr
01972     );
01973 }
01974 
01975 
01976 bool Foam::polyMeshGeometry::checkFaceWeights
01977 (
01978     const bool report,
01979     const scalar warnWeight,
01980     const labelList& checkFaces,
01981     const List<labelPair>& baffles,
01982     labelHashSet* setPtr
01983 ) const
01984 {
01985     return checkFaceWeights
01986     (
01987         report,
01988         warnWeight,
01989         mesh_,
01990         cellCentres_,
01991         faceCentres_,
01992         faceAreas_,
01993         checkFaces,
01994         baffles,
01995         setPtr
01996     );
01997 }
01998 
01999 
02000 bool Foam::polyMeshGeometry::checkVolRatio
02001 (
02002     const bool report,
02003     const scalar warnRatio,
02004     const labelList& checkFaces,
02005     const List<labelPair>& baffles,
02006     labelHashSet* setPtr
02007 ) const
02008 {
02009     return checkVolRatio
02010     (
02011         report,
02012         warnRatio,
02013         mesh_,
02014         cellVolumes_,
02015         checkFaces,
02016         baffles,
02017         setPtr
02018     );
02019 }
02020 
02021 
02022 bool Foam::polyMeshGeometry::checkFaceAngles
02023 (
02024     const bool report,
02025     const scalar maxDeg,
02026     const pointField& p,
02027     const labelList& checkFaces,
02028     labelHashSet* setPtr
02029 ) const
02030 {
02031     return checkFaceAngles
02032     (
02033         report,
02034         maxDeg,
02035         mesh_,
02036         faceAreas_,
02037         p,
02038         checkFaces,
02039         setPtr
02040     );
02041 }
02042 
02043 
02044 bool Foam::polyMeshGeometry::checkFaceTwist
02045 (
02046     const bool report,
02047     const scalar minTwist,
02048     const pointField& p,
02049     const labelList& checkFaces,
02050     labelHashSet* setPtr
02051 ) const
02052 {
02053     return checkFaceTwist
02054     (
02055         report,
02056         minTwist,
02057         mesh_,
02058         cellCentres_,
02059         faceAreas_,
02060         faceCentres_,
02061         p,
02062         checkFaces,
02063         setPtr
02064     );
02065 }
02066 
02067 
02068 bool Foam::polyMeshGeometry::checkTriangleTwist
02069 (
02070     const bool report,
02071     const scalar minTwist,
02072     const pointField& p,
02073     const labelList& checkFaces,
02074     labelHashSet* setPtr
02075 ) const
02076 {
02077     return checkTriangleTwist
02078     (
02079         report,
02080         minTwist,
02081         mesh_,
02082         faceAreas_,
02083         faceCentres_,
02084         p,
02085         checkFaces,
02086         setPtr
02087     );
02088 }
02089 
02090 
02091 bool Foam::polyMeshGeometry::checkFaceArea
02092 (
02093     const bool report,
02094     const scalar minArea,
02095     const labelList& checkFaces,
02096     labelHashSet* setPtr
02097 ) const
02098 {
02099     return checkFaceArea
02100     (
02101         report,
02102         minArea,
02103         mesh_,
02104         faceAreas_,
02105         checkFaces,
02106         setPtr
02107     );
02108 }
02109 
02110 
02111 bool Foam::polyMeshGeometry::checkCellDeterminant
02112 (
02113     const bool report,
02114     const scalar warnDet,
02115     const labelList& checkFaces,
02116     const labelList& affectedCells,
02117     labelHashSet* setPtr
02118 ) const
02119 {
02120     return checkCellDeterminant
02121     (
02122         report,
02123         warnDet,
02124         mesh_,
02125         faceAreas_,
02126         checkFaces,
02127         affectedCells,
02128         setPtr
02129     );
02130 }
02131 
02132 
02133 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines