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

primitiveMeshGeometry.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 "primitiveMeshGeometry.H"
00027 #include <OpenFOAM/pyramidPointFaceRef.H>
00028 
00029 namespace Foam
00030 {
00031 
00032 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00033 
00034 defineTypeNameAndDebug(primitiveMeshGeometry, 0);
00035 
00036 }
00037 
00038 
00039 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00040 
00041 void Foam::primitiveMeshGeometry::updateFaceCentresAndAreas
00042 (
00043     const pointField& p,
00044     const labelList& changedFaces
00045 )
00046 {
00047     const faceList& fs = mesh_.faces();
00048 
00049     forAll(changedFaces, i)
00050     {
00051         label facei = changedFaces[i];
00052 
00053         const labelList& f = fs[facei];
00054         label nPoints = f.size();
00055 
00056         // If the face is a triangle, do a direct calculation for efficiency
00057         // and to avoid round-off error-related problems
00058         if (nPoints == 3)
00059         {
00060             faceCentres_[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
00061             faceAreas_[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
00062         }
00063         else
00064         {
00065             vector sumN = vector::zero;
00066             scalar sumA = 0.0;
00067             vector sumAc = vector::zero;
00068 
00069             point fCentre = p[f[0]];
00070             for (label pi = 1; pi < nPoints; pi++)
00071             {
00072                 fCentre += p[f[pi]];
00073             }
00074 
00075             fCentre /= nPoints;
00076 
00077             for (label pi = 0; pi < nPoints; pi++)
00078             {
00079                 const point& nextPoint = p[f[(pi + 1) % nPoints]];
00080 
00081                 vector c = p[f[pi]] + nextPoint + fCentre;
00082                 vector n = (nextPoint - p[f[pi]])^(fCentre - p[f[pi]]);
00083                 scalar a = mag(n);
00084 
00085                 sumN += n;
00086                 sumA += a;
00087                 sumAc += a*c;
00088             }
00089 
00090             faceCentres_[facei] = (1.0/3.0)*sumAc/(sumA + VSMALL);
00091             faceAreas_[facei] = 0.5*sumN;
00092         }
00093     }
00094 }
00095 
00096 
00097 void Foam::primitiveMeshGeometry::updateCellCentresAndVols
00098 (
00099     const labelList& changedCells,
00100     const labelList& changedFaces
00101 )
00102 {
00103     // Clear the fields for accumulation
00104     UIndirectList<vector>(cellCentres_, changedCells) = vector::zero;
00105     UIndirectList<scalar>(cellVolumes_, changedCells) = 0.0;
00106 
00107     const labelList& own = mesh_.faceOwner();
00108     const labelList& nei = mesh_.faceNeighbour();
00109 
00110     // first estimate the approximate cell centre as the average of face centres
00111 
00112     vectorField cEst(mesh_.nCells());
00113     UIndirectList<vector>(cEst, changedCells) = vector::zero;
00114     scalarField nCellFaces(mesh_.nCells());
00115     UIndirectList<scalar>(nCellFaces, changedCells) = 0.0;
00116 
00117     forAll(changedFaces, i)
00118     {
00119         label faceI = changedFaces[i];
00120         cEst[own[faceI]] += faceCentres_[faceI];
00121         nCellFaces[own[faceI]] += 1;
00122 
00123         if (mesh_.isInternalFace(faceI))
00124         {
00125             cEst[nei[faceI]] += faceCentres_[faceI];
00126             nCellFaces[nei[faceI]] += 1;
00127         }
00128     }
00129 
00130     forAll(changedCells, i)
00131     {
00132         label cellI = changedCells[i];
00133         cEst[cellI] /= nCellFaces[cellI];
00134     }
00135 
00136     forAll(changedFaces, i)
00137     {
00138         label faceI = changedFaces[i];
00139 
00140         // Calculate 3*face-pyramid volume
00141         scalar pyr3Vol = max
00142         (
00143             faceAreas_[faceI] & (faceCentres_[faceI] - cEst[own[faceI]]),
00144             VSMALL
00145         );
00146 
00147         // Calculate face-pyramid centre
00148         vector pc = (3.0/4.0)*faceCentres_[faceI] + (1.0/4.0)*cEst[own[faceI]];
00149 
00150         // Accumulate volume-weighted face-pyramid centre
00151         cellCentres_[own[faceI]] += pyr3Vol*pc;
00152 
00153         // Accumulate face-pyramid volume
00154         cellVolumes_[own[faceI]] += pyr3Vol;
00155 
00156         if (mesh_.isInternalFace(faceI))
00157         {
00158             // Calculate 3*face-pyramid volume
00159             scalar pyr3Vol = max
00160             (
00161                 faceAreas_[faceI] & (cEst[nei[faceI]] - faceCentres_[faceI]),
00162                 VSMALL
00163             );
00164 
00165             // Calculate face-pyramid centre
00166             vector pc =
00167                 (3.0/4.0)*faceCentres_[faceI]
00168               + (1.0/4.0)*cEst[nei[faceI]];
00169 
00170             // Accumulate volume-weighted face-pyramid centre
00171             cellCentres_[nei[faceI]] += pyr3Vol*pc;
00172 
00173             // Accumulate face-pyramid volume
00174             cellVolumes_[nei[faceI]] += pyr3Vol;
00175         }
00176     }
00177 
00178     forAll(changedCells, i)
00179     {
00180         label cellI = changedCells[i];
00181 
00182         cellCentres_[cellI] /= cellVolumes_[cellI];
00183         cellVolumes_[cellI] *= (1.0/3.0);
00184     }
00185 }
00186 
00187 
00188 Foam::labelList Foam::primitiveMeshGeometry::affectedCells
00189 (
00190     const labelList& changedFaces
00191 ) const
00192 {
00193     const labelList& own = mesh_.faceOwner();
00194     const labelList& nei = mesh_.faceNeighbour();
00195 
00196     labelHashSet affectedCells(2*changedFaces.size());
00197 
00198     forAll(changedFaces, i)
00199     {
00200         label faceI = changedFaces[i];
00201 
00202         affectedCells.insert(own[faceI]);
00203 
00204         if (mesh_.isInternalFace(faceI))
00205         {
00206             affectedCells.insert(nei[faceI]);
00207         }
00208     }
00209     return affectedCells.toc();
00210 }
00211 
00212 
00213 
00214 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00215 
00216 // Construct from components
00217 Foam::primitiveMeshGeometry::primitiveMeshGeometry
00218 (
00219     const primitiveMesh& mesh
00220 )
00221 :
00222     mesh_(mesh)
00223 {
00224     correct();
00225 }
00226 
00227 
00228 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00229 
00230 
00231 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00232 
00233 //- Take over properties from mesh
00234 void Foam::primitiveMeshGeometry::correct()
00235 {
00236     faceAreas_ = mesh_.faceAreas();
00237     faceCentres_ = mesh_.faceCentres();
00238     cellCentres_ = mesh_.cellCentres();
00239     cellVolumes_ = mesh_.cellVolumes();
00240 }
00241 
00242 
00243 //- Recalculate on selected faces
00244 void Foam::primitiveMeshGeometry::correct
00245 (
00246     const pointField& p,
00247     const labelList& changedFaces
00248 )
00249 {
00250     // Update face quantities
00251     updateFaceCentresAndAreas(p, changedFaces);
00252     // Update cell quantities from face quantities
00253     updateCellCentresAndVols(affectedCells(changedFaces), changedFaces);
00254 }
00255 
00256 
00257 bool Foam::primitiveMeshGeometry::checkFaceDotProduct
00258 (
00259     const bool report,
00260     const scalar orthWarn,
00261     const primitiveMesh& mesh,
00262     const vectorField& cellCentres,
00263     const vectorField& faceAreas,
00264     const labelList& checkFaces,
00265     labelHashSet* setPtr
00266 )
00267 {
00268     // for all internal faces check theat the d dot S product is positive
00269 
00270     const labelList& own = mesh.faceOwner();
00271     const labelList& nei = mesh.faceNeighbour();
00272 
00273     // Severe nonorthogonality threshold
00274     const scalar severeNonorthogonalityThreshold =
00275         ::cos(orthWarn/180.0*mathematicalConstant::pi);
00276 
00277     scalar minDDotS = GREAT;
00278 
00279     scalar sumDDotS = 0;
00280 
00281     label severeNonOrth = 0;
00282 
00283     label errorNonOrth = 0;
00284 
00285     forAll(checkFaces, i)
00286     {
00287         label faceI = checkFaces[i];
00288 
00289         if (mesh.isInternalFace(faceI))
00290         {
00291             vector d = cellCentres[nei[faceI]] - cellCentres[own[faceI]];
00292             const vector& s = faceAreas[faceI];
00293 
00294             scalar dDotS = (d & s)/(mag(d)*mag(s) + VSMALL);
00295 
00296             if (dDotS < severeNonorthogonalityThreshold)
00297             {
00298                 if (dDotS > SMALL)
00299                 {
00300                     if (report)
00301                     {
00302                         // Severe non-orthogonality but mesh still OK
00303                         Pout<< "Severe non-orthogonality for face " << faceI
00304                             << " between cells " << own[faceI]
00305                             << " and " << nei[faceI]
00306                             << ": Angle = "
00307                             << ::acos(dDotS)/mathematicalConstant::pi*180.0
00308                             << " deg." << endl;
00309                     }
00310 
00311                     if (setPtr)
00312                     {
00313                         setPtr->insert(faceI);
00314                     }
00315 
00316                     severeNonOrth++;
00317                 }
00318                 else
00319                 {
00320                     // Non-orthogonality greater than 90 deg
00321                     if (report)
00322                     {
00323                         WarningIn
00324                         (
00325                             "primitiveMeshGeometry::checkFaceDotProduct"
00326                             "(const bool, const scalar, const labelList&"
00327                             ", labelHashSet*)"
00328                         )   << "Severe non-orthogonality detected for face "
00329                             << faceI
00330                             << " between cells " << own[faceI] << " and "
00331                             << nei[faceI]
00332                             << ": Angle = "
00333                             << ::acos(dDotS)/mathematicalConstant::pi*180.0
00334                             << " deg." << endl;
00335                     }
00336 
00337                     errorNonOrth++;
00338 
00339                     if (setPtr)
00340                     {
00341                         setPtr->insert(faceI);
00342                     }
00343                 }
00344             }
00345 
00346             if (dDotS < minDDotS)
00347             {
00348                 minDDotS = dDotS;
00349             }
00350 
00351             sumDDotS += dDotS;
00352         }
00353     }
00354 
00355     reduce(minDDotS, minOp<scalar>());
00356     reduce(sumDDotS, sumOp<scalar>());
00357     reduce(severeNonOrth, sumOp<label>());
00358     reduce(errorNonOrth, sumOp<label>());
00359 
00360     label neiSize = nei.size();
00361     reduce(neiSize, sumOp<label>());
00362 
00363     // Only report if there are some internal faces
00364     if (neiSize > 0)
00365     {
00366         if (report && minDDotS < severeNonorthogonalityThreshold)
00367         {
00368             Info<< "Number of non-orthogonality errors: " << errorNonOrth
00369                 << ". Number of severely non-orthogonal faces: "
00370                 << severeNonOrth  << "." << endl;
00371         }
00372     }
00373 
00374     if (report)
00375     {
00376         if (neiSize > 0)
00377         {
00378             Info<< "Mesh non-orthogonality Max: "
00379                 << ::acos(minDDotS)/mathematicalConstant::pi*180.0
00380                 << " average: " <<
00381                    ::acos(sumDDotS/neiSize)/mathematicalConstant::pi*180.0
00382                 << endl;
00383         }
00384     }
00385 
00386     if (errorNonOrth > 0)
00387     {
00388         if (report)
00389         {
00390             SeriousErrorIn
00391             (
00392                 "primitiveMeshGeometry::checkFaceDotProduct"
00393                 "(const bool, const scalar, const labelList&, labelHashSet*)"
00394             )   << "Error in non-orthogonality detected" << endl;
00395         }
00396 
00397         return true;
00398     }
00399     else
00400     {
00401         if (report)
00402         {
00403             Info<< "Non-orthogonality check OK.\n" << endl;
00404         }
00405 
00406         return false;
00407     }
00408 }
00409 
00410 
00411 bool Foam::primitiveMeshGeometry::checkFacePyramids
00412 (
00413     const bool report,
00414     const scalar minPyrVol,
00415     const primitiveMesh& mesh,
00416     const vectorField& cellCentres,
00417     const pointField& p,
00418     const labelList& checkFaces,
00419     labelHashSet* setPtr
00420 )
00421 {
00422     // check whether face area vector points to the cell with higher label
00423     const labelList& own = mesh.faceOwner();
00424     const labelList& nei = mesh.faceNeighbour();
00425 
00426     const faceList& f = mesh.faces();
00427 
00428     label nErrorPyrs = 0;
00429 
00430     forAll(checkFaces, i)
00431     {
00432         label faceI = checkFaces[i];
00433 
00434         // Create the owner pyramid - it will have negative volume
00435         scalar pyrVol = pyramidPointFaceRef
00436         (
00437             f[faceI],
00438             cellCentres[own[faceI]]
00439         ).mag(p);
00440 
00441         if (pyrVol > -minPyrVol)
00442         {
00443             if (report)
00444             {
00445                 Pout<< "bool primitiveMeshGeometry::checkFacePyramids("
00446                     << "const bool, const scalar, const pointField&"
00447                     << ", const labelList&, labelHashSet*): "
00448                     << "face " << faceI << " points the wrong way. " << endl
00449                     << "Pyramid volume: " << -pyrVol
00450                     << " Face " << f[faceI] << " area: " << f[faceI].mag(p)
00451                     << " Owner cell: " << own[faceI] << endl
00452                     << "Owner cell vertex labels: "
00453                     << mesh.cells()[own[faceI]].labels(f)
00454                     << endl;
00455             }
00456 
00457 
00458             if (setPtr)
00459             {
00460                 setPtr->insert(faceI);
00461             }
00462 
00463             nErrorPyrs++;
00464         }
00465 
00466         if (mesh.isInternalFace(faceI))
00467         {
00468             // Create the neighbour pyramid - it will have positive volume
00469             scalar pyrVol =
00470                 pyramidPointFaceRef(f[faceI], cellCentres[nei[faceI]]).mag(p);
00471 
00472             if (pyrVol < minPyrVol)
00473             {
00474                 if (report)
00475                 {
00476                     Pout<< "bool primitiveMeshGeometry::checkFacePyramids("
00477                         << "const bool, const scalar, const pointField&"
00478                         << ", const labelList&, labelHashSet*): "
00479                         << "face " << faceI << " points the wrong way. " << endl
00480                         << "Pyramid volume: " << -pyrVol
00481                         << " Face " << f[faceI] << " area: " << f[faceI].mag(p)
00482                         << " Neighbour cell: " << nei[faceI] << endl
00483                         << "Neighbour cell vertex labels: "
00484                         << mesh.cells()[nei[faceI]].labels(f)
00485                         << endl;
00486                 }
00487 
00488                 if (setPtr)
00489                 {
00490                     setPtr->insert(faceI);
00491                 }
00492 
00493                 nErrorPyrs++;
00494             }
00495         }
00496     }
00497 
00498     reduce(nErrorPyrs, sumOp<label>());
00499 
00500     if (nErrorPyrs > 0)
00501     {
00502         if (report)
00503         {
00504             SeriousErrorIn
00505             (
00506                 "primitiveMeshGeometry::checkFacePyramids("
00507                 "const bool, const scalar, const pointField&"
00508                 ", const labelList&, labelHashSet*)"
00509             )   << "Error in face pyramids: faces pointing the wrong way!"
00510                 << endl;
00511         }
00512 
00513         return true;
00514     }
00515     else
00516     {
00517         if (report)
00518         {
00519             Info<< "Face pyramids OK.\n" << endl;
00520         }
00521 
00522         return false;
00523     }
00524 }
00525 
00526 
00527 bool Foam::primitiveMeshGeometry::checkFaceSkewness
00528 (
00529     const bool report,
00530     const scalar internalSkew,
00531     const scalar boundarySkew,
00532     const primitiveMesh& mesh,
00533     const vectorField& cellCentres,
00534     const vectorField& faceCentres,
00535     const vectorField& faceAreas,
00536     const labelList& checkFaces,
00537     labelHashSet* setPtr
00538 )
00539 {
00540     // Warn if the skew correction vector is more than skew times
00541     // larger than the face area vector
00542 
00543     const labelList& own = mesh.faceOwner();
00544     const labelList& nei = mesh.faceNeighbour();
00545 
00546     scalar maxSkew = 0;
00547 
00548     label nWarnSkew = 0;
00549 
00550     forAll(checkFaces, i)
00551     {
00552         label faceI = checkFaces[i];
00553 
00554         if (mesh.isInternalFace(faceI))
00555         {
00556             scalar dOwn = mag(faceCentres[faceI] - cellCentres[own[faceI]]);
00557             scalar dNei = mag(faceCentres[faceI] - cellCentres[nei[faceI]]);
00558 
00559             point faceIntersection =
00560                 cellCentres[own[faceI]]*dNei/(dOwn+dNei)
00561               + cellCentres[nei[faceI]]*dOwn/(dOwn+dNei);
00562 
00563             scalar skewness = 
00564                 mag(faceCentres[faceI] - faceIntersection)
00565               / (
00566                     mag(cellCentres[nei[faceI]]-cellCentres[own[faceI]])
00567                   + VSMALL
00568                 );
00569 
00570             // Check if the skewness vector is greater than the PN vector.
00571             // This does not cause trouble but is a good indication of a poor
00572             // mesh.
00573             if (skewness > internalSkew)
00574             {
00575                 if (report)
00576                 {
00577                     Pout<< "Severe skewness for face " << faceI
00578                         << " skewness = " << skewness << endl;
00579                 }
00580 
00581                 if (setPtr)
00582                 {
00583                     setPtr->insert(faceI);
00584                 }
00585 
00586                 nWarnSkew++;
00587             }
00588 
00589             if (skewness > maxSkew)
00590             {
00591                 maxSkew = skewness;
00592             }
00593         }
00594         else
00595         {
00596             // Boundary faces: consider them to have only skewness error.
00597             // (i.e. treat as if mirror cell on other side)
00598 
00599             vector faceNormal = faceAreas[faceI];
00600             faceNormal /= mag(faceNormal) + VSMALL;
00601 
00602             vector dOwn = faceCentres[faceI] - cellCentres[own[faceI]];
00603 
00604             vector dWall = faceNormal*(faceNormal & dOwn);
00605 
00606             point faceIntersection = cellCentres[own[faceI]] + dWall;
00607 
00608             scalar skewness =
00609                 mag(faceCentres[faceI] - faceIntersection)
00610                 /(2*mag(dWall) + VSMALL);
00611 
00612             // Check if the skewness vector is greater than the PN vector.
00613             // This does not cause trouble but is a good indication of a poor
00614             // mesh.
00615             if (skewness > boundarySkew)
00616             {
00617                 if (report)
00618                 {
00619                     Pout<< "Severe skewness for boundary face " << faceI
00620                         << " skewness = " << skewness << endl;
00621                 }
00622 
00623                 if (setPtr)
00624                 {
00625                     setPtr->insert(faceI);
00626                 }
00627 
00628                 nWarnSkew++;
00629             }
00630 
00631             if (skewness > maxSkew)
00632             {
00633                 maxSkew = skewness;
00634             }
00635         }
00636     }
00637 
00638     reduce(maxSkew, maxOp<scalar>());
00639     reduce(nWarnSkew, sumOp<label>());
00640 
00641     if (nWarnSkew > 0)
00642     {
00643         if (report)
00644         {
00645             WarningIn
00646             (
00647                 "primitiveMeshGeometry::checkFaceSkewness"
00648                 "(const bool, const scalar, const labelList&, labelHashSet*)"
00649             )   << "Large face skewness detected.  Max skewness = "
00650                 << 100*maxSkew
00651                 << " percent.\nThis may impair the quality of the result." << nl
00652                 << nWarnSkew << " highly skew faces detected."
00653                 << endl;
00654         }
00655 
00656         return true;
00657     }
00658     else
00659     {
00660         if (report)
00661         {
00662             Info<< "Max skewness = " << 100*maxSkew
00663                 << " percent.  Face skewness OK.\n" << endl;
00664         }
00665 
00666         return false;
00667     }
00668 }
00669 
00670 
00671 bool Foam::primitiveMeshGeometry::checkFaceWeights
00672 (
00673     const bool report,
00674     const scalar warnWeight,
00675     const primitiveMesh& mesh,
00676     const vectorField& cellCentres,
00677     const vectorField& faceCentres,
00678     const vectorField& faceAreas,
00679     const labelList& checkFaces,
00680     labelHashSet* setPtr
00681 )
00682 {
00683     // Warn if the delta factor (0..1) is too large.
00684 
00685     const labelList& own = mesh.faceOwner();
00686     const labelList& nei = mesh.faceNeighbour();
00687 
00688     scalar minWeight = GREAT;
00689 
00690     label nWarnWeight = 0;
00691 
00692     forAll(checkFaces, i)
00693     {
00694         label faceI = checkFaces[i];
00695 
00696         if (mesh.isInternalFace(faceI))
00697         {
00698             const point& fc = faceCentres[faceI];
00699 
00700             scalar dOwn = mag(faceAreas[faceI] & (fc-cellCentres[own[faceI]]));
00701             scalar dNei = mag(faceAreas[faceI] & (cellCentres[nei[faceI]]-fc));
00702 
00703             scalar weight = min(dNei,dOwn)/(dNei+dOwn);
00704 
00705             if (weight < warnWeight)
00706             {
00707                 if (report)
00708                 {
00709                     Pout<< "Small weighting factor for face " << faceI
00710                         << " weight = " << weight << endl;
00711                 }
00712 
00713                 if (setPtr)
00714                 {
00715                     setPtr->insert(faceI);
00716                 }
00717 
00718                 nWarnWeight++;
00719             }
00720 
00721             minWeight = min(minWeight, weight);
00722         }
00723     }
00724 
00725     reduce(minWeight, minOp<scalar>());
00726     reduce(nWarnWeight, sumOp<label>());
00727 
00728     if (minWeight < warnWeight)
00729     {
00730         if (report)
00731         {
00732             WarningIn
00733             (
00734                 "primitiveMeshGeometry::checkFaceWeights"
00735                 "(const bool, const scalar, const labelList&, labelHashSet*)"
00736             )   << "Small interpolation weight detected.  Min weight = "
00737                 << minWeight << '.' << nl
00738                 << nWarnWeight << " faces with small weights detected."
00739                 << endl;
00740         }
00741 
00742         return true;
00743     }
00744     else
00745     {
00746         if (report)
00747         {
00748             Info<< "Min weight = " << minWeight
00749                 << " percent.  Weights OK.\n" << endl;
00750         }
00751 
00752         return false;
00753     }
00754 }
00755 
00756 
00757 // Check convexity of angles in a face. Allow a slight non-convexity.
00758 // E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity)
00759 // (if truly concave and points not visible from face centre the face-pyramid
00760 //  check in checkMesh will fail)
00761 bool Foam::primitiveMeshGeometry::checkFaceAngles
00762 (
00763     const bool report,
00764     const scalar maxDeg,
00765     const primitiveMesh& mesh,
00766     const vectorField& faceAreas,
00767     const pointField& p,
00768     const labelList& checkFaces,
00769     labelHashSet* setPtr
00770 )
00771 {
00772     if (maxDeg < -SMALL || maxDeg > 180+SMALL)
00773     {
00774         FatalErrorIn
00775         (
00776             "primitiveMeshGeometry::checkFaceAngles"
00777             "(const bool, const scalar, const pointField&, const labelList&"
00778             ", labelHashSet*)"
00779         )   << "maxDeg should be [0..180] but is now " << maxDeg
00780             << abort(FatalError);
00781     }
00782 
00783     const scalar maxSin = Foam::sin(maxDeg/180.0*mathematicalConstant::pi);
00784 
00785     const faceList& fcs = mesh.faces();
00786 
00787     scalar maxEdgeSin = 0.0;
00788 
00789     label nConcave = 0;
00790 
00791     label errorFaceI = -1;
00792 
00793     forAll(checkFaces, i)
00794     {
00795         label faceI = checkFaces[i];
00796 
00797         const face& f = fcs[faceI];
00798 
00799         vector faceNormal = faceAreas[faceI];
00800         faceNormal /= mag(faceNormal) + VSMALL;
00801 
00802         // Get edge from f[0] to f[size-1];
00803         vector ePrev(p[f[0]] - p[f[f.size()-1]]);
00804         scalar magEPrev = mag(ePrev);
00805         ePrev /= magEPrev + VSMALL;
00806 
00807         forAll(f, fp0)
00808         {
00809             // Get vertex after fp
00810             label fp1 = f.fcIndex(fp0);
00811 
00812             // Normalized vector between two consecutive points
00813             vector e10(p[f[fp1]] - p[f[fp0]]);
00814             scalar magE10 = mag(e10);
00815             e10 /= magE10 + VSMALL;
00816 
00817             if (magEPrev > SMALL && magE10 > SMALL)
00818             {
00819                 vector edgeNormal = ePrev ^ e10;
00820                 scalar magEdgeNormal = mag(edgeNormal);
00821 
00822                 if (magEdgeNormal < maxSin)
00823                 {
00824                     // Edges (almost) aligned -> face is ok.
00825                 }
00826                 else
00827                 {
00828                     // Check normal
00829                     edgeNormal /= magEdgeNormal;
00830 
00831                     if ((edgeNormal & faceNormal) < SMALL)
00832                     {
00833                         if (faceI != errorFaceI)
00834                         {
00835                             // Count only one error per face.
00836                             errorFaceI = faceI;
00837                             nConcave++;
00838                         }
00839 
00840                         if (setPtr)
00841                         {
00842                             setPtr->insert(faceI);
00843                         }
00844 
00845                         maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
00846                     }
00847                 }
00848             }
00849 
00850             ePrev = e10;
00851             magEPrev = magE10;
00852         }
00853     }
00854 
00855     reduce(nConcave, sumOp<label>());
00856     reduce(maxEdgeSin, maxOp<scalar>());
00857 
00858     if (report)
00859     {
00860         if (maxEdgeSin > SMALL)
00861         {
00862             scalar maxConcaveDegr =
00863                 Foam::asin(Foam::min(1.0, maxEdgeSin))
00864              * 180.0/mathematicalConstant::pi;
00865 
00866             Info<< "There are " << nConcave
00867                 << " faces with concave angles between consecutive"
00868                 << " edges. Max concave angle = " << maxConcaveDegr
00869                 << " degrees.\n" << endl;
00870         }
00871         else
00872         {
00873             Info<< "All angles in faces are convex or less than "  << maxDeg
00874                 << " degrees concave.\n" << endl;
00875         }
00876     }
00877 
00878     if (nConcave > 0)
00879     {
00880         if (report)
00881         {
00882             WarningIn
00883             (
00884                 "primitiveMeshGeometry::checkFaceAngles"
00885                 "(const bool, const scalar,  const pointField&"
00886                 ", const labelList&, labelHashSet*)"
00887             )   << nConcave  << " face points with severe concave angle (> "
00888                 << maxDeg << " deg) found.\n"
00889                 << endl;
00890         }
00891 
00892         return true;
00893     }
00894     else
00895     {
00896         return false;
00897     }
00898 }
00899 
00900 
00904 //bool Foam::primitiveMeshGeometry::checkFaceFlatness
00905 //(
00906 //    const bool report,
00907 //    const scalar warnFlatness,
00908 //    const primitiveMesh& mesh,
00909 //    const vectorField& faceAreas,
00910 //    const vectorField& faceCentres,
00911 //    const pointField& p,
00912 //    const labelList& checkFaces,
00913 //    labelHashSet* setPtr
00914 //)
00915 //{
00916 //    if (warnFlatness < 0 || warnFlatness > 1)
00917 //    {
00918 //        FatalErrorIn
00919 //        (
00920 //            "primitiveMeshGeometry::checkFaceFlatness"
00921 //            "(const bool, const scalar, const pointField&"
00922 //            ", const labelList&, labelHashSet*)"
00923 //        )   << "warnFlatness should be [0..1] but is now " << warnFlatness
00924 //            << abort(FatalError);
00925 //    }
00926 //
00927 //
00928 //    const faceList& fcs = mesh.faces();
00929 //
00930 //    // Areas are calculated as the sum of areas. (see
00931 //    // primitiveMeshFaceCentresAndAreas.C)
00932 //
00933 //    label nWarped = 0;
00934 //
00935 //    scalar minFlatness = GREAT;
00936 //    scalar sumFlatness = 0;
00937 //    label nSummed = 0;
00938 //
00939 //    forAll(checkFaces, i)
00940 //    {
00941 //        label faceI = checkFaces[i];
00942 //
00943 //        const face& f = fcs[faceI];
00944 //
00945 //        scalar magArea = mag(faceAreas[faceI]);
00946 //
00947 //        if (f.size() > 3 && magArea > VSMALL)
00948 //        {
00949 //            const point& fc = faceCentres[faceI];
00950 //
00951 //            // Calculate the sum of magnitude of areas and compare to magnitude
00952 //            // of sum of areas.
00953 //
00954 //            scalar sumA = 0.0;
00955 //
00956 //            forAll(f, fp)
00957 //            {
00958 //                const point& thisPoint = p[f[fp]];
00959 //                const point& nextPoint = p[f.nextLabel(fp)];
00960 //
00961 //                // Triangle around fc.
00962 //                vector n = 0.5*((nextPoint - thisPoint)^(fc - thisPoint));
00963 //                sumA += mag(n);
00964 //            }
00965 //
00966 //            scalar flatness = magArea / (sumA+VSMALL);
00967 //
00968 //            sumFlatness += flatness;
00969 //            nSummed++;
00970 //
00971 //            minFlatness = min(minFlatness, flatness);
00972 //
00973 //            if (flatness < warnFlatness)
00974 //            {
00975 //                nWarped++;
00976 //
00977 //                if (setPtr)
00978 //                {
00979 //                    setPtr->insert(faceI);
00980 //                }
00981 //            }
00982 //        }
00983 //    }
00984 //
00985 //
00986 //    reduce(nWarped, sumOp<label>());
00987 //    reduce(minFlatness, minOp<scalar>());
00988 //
00989 //    reduce(nSummed, sumOp<label>());
00990 //    reduce(sumFlatness, sumOp<scalar>());
00991 //
00992 //    if (report)
00993 //    {
00994 //        if (nSummed > 0)
00995 //        {
00996 //            Info<< "Face flatness (1 = flat, 0 = butterfly) : average = "
00997 //                << sumFlatness / nSummed << "  min = " << minFlatness << endl;
00998 //        }
00999 //
01000 //        if (nWarped> 0)
01001 //        {
01002 //            Info<< "There are " << nWarped
01003 //                << " faces with ratio between projected and actual area < "
01004 //                << warnFlatness
01005 //                << ".\nMinimum ratio (minimum flatness, maximum warpage) = "
01006 //                << minFlatness << nl << endl;
01007 //        }
01008 //        else
01009 //        {
01010 //            Info<< "All faces are flat in that the ratio between projected"
01011 //                << " and actual area is > " << warnFlatness << nl << endl;
01012 //        }
01013 //    }
01014 //
01015 //    if (nWarped > 0)
01016 //    {
01017 //        if (report)
01018 //        {
01019 //            WarningIn
01020 //            (
01021 //                "primitiveMeshGeometry::checkFaceFlatness"
01022 //                "(const bool, const scalar, const pointField&"
01023 //                ", const labelList&, labelHashSet*)"
01024 //            )   << nWarped  << " faces with severe warpage (flatness < "
01025 //                << warnFlatness << ") found.\n"
01026 //                << endl;
01027 //        }
01028 //
01029 //        return true;
01030 //    }
01031 //    else
01032 //    {
01033 //        return false;
01034 //    }
01035 //}
01036 
01037 
01038 // Check twist of faces. Is calculated as the difference between areas of
01039 // individual triangles and the overall area of the face (which ifself is
01040 // is the average of the areas of the individual triangles).
01041 bool Foam::primitiveMeshGeometry::checkFaceTwist
01042 (
01043     const bool report,
01044     const scalar minTwist,
01045     const primitiveMesh& mesh,
01046     const vectorField& faceAreas,
01047     const vectorField& faceCentres,
01048     const pointField& p,
01049     const labelList& checkFaces,
01050     labelHashSet* setPtr
01051 )
01052 {
01053     if (minTwist < -1-SMALL || minTwist > 1+SMALL)
01054     {
01055         FatalErrorIn
01056         (
01057             "primitiveMeshGeometry::checkFaceTwist"
01058             "(const bool, const scalar, const primitiveMesh&, const pointField&"
01059             ", const labelList&, labelHashSet*)"
01060         )   << "minTwist should be [-1..1] but is now " << minTwist
01061             << abort(FatalError);
01062     }
01063 
01064 
01065     const faceList& fcs = mesh.faces();
01066 
01067     // Areas are calculated as the sum of areas. (see
01068     // primitiveMeshFaceCentresAndAreas.C)
01069 
01070     label nWarped = 0;
01071 
01072     forAll(checkFaces, i)
01073     {
01074         label faceI = checkFaces[i];
01075 
01076         const face& f = fcs[faceI];
01077 
01078         scalar magArea = mag(faceAreas[faceI]);
01079 
01080         if (f.size() > 3 && magArea > VSMALL)
01081         {
01082             const vector nf = faceAreas[faceI] / magArea;
01083 
01084             const point& fc = faceCentres[faceI];
01085 
01086             forAll(f, fpI)
01087             {
01088                 vector triArea
01089                 (
01090                     triPointRef
01091                     (
01092                         p[f[fpI]],
01093                         p[f.nextLabel(fpI)],
01094                         fc
01095                     ).normal()
01096                 );
01097 
01098                 scalar magTri = mag(triArea);
01099 
01100                 if (magTri > VSMALL && ((nf & triArea/magTri) < minTwist))
01101                 {
01102                     nWarped++;
01103 
01104                     if (setPtr)
01105                     {
01106                         setPtr->insert(faceI);
01107                     }
01108                 }
01109             }
01110         }
01111     }
01112 
01113 
01114     reduce(nWarped, sumOp<label>());
01115 
01116     if (report)
01117     {
01118         if (nWarped> 0)
01119         {
01120             Info<< "There are " << nWarped
01121                 << " faces with cosine of the angle"
01122                 << " between triangle normal and face normal less than "
01123                 << minTwist << nl << endl;
01124         }
01125         else
01126         {
01127             Info<< "All faces are flat in that the cosine of the angle"
01128                 << " between triangle normal and face normal less than "
01129                 << minTwist << nl << endl;
01130         }
01131     }
01132 
01133     if (nWarped > 0)
01134     {
01135         if (report)
01136         {
01137             WarningIn
01138             (
01139                 "primitiveMeshGeometry::checkFaceTwist"
01140                 "(const bool, const scalar, const primitiveMesh&"
01141                 ", const pointField&, const labelList&, labelHashSet*)"
01142             )   << nWarped  << " faces with severe warpage "
01143                 << "(cosine of the angle between triangle normal and face normal"
01144                 << " < " << minTwist << ") found.\n"
01145                 << endl;
01146         }
01147 
01148         return true;
01149     }
01150     else
01151     {
01152         return false;
01153     }
01154 }
01155 
01156 
01157 bool Foam::primitiveMeshGeometry::checkFaceArea
01158 (
01159     const bool report,
01160     const scalar minArea,
01161     const primitiveMesh& mesh,
01162     const vectorField& faceAreas,
01163     const labelList& checkFaces,
01164     labelHashSet* setPtr
01165 )
01166 {
01167     label nZeroArea = 0;
01168 
01169     forAll(checkFaces, i)
01170     {
01171         label faceI = checkFaces[i];
01172 
01173         if (mag(faceAreas[faceI]) < minArea)
01174         {
01175             if (setPtr)
01176             {
01177                 setPtr->insert(faceI);
01178             }
01179             nZeroArea++;
01180         }
01181     }
01182 
01183 
01184     reduce(nZeroArea, sumOp<label>());
01185 
01186     if (report)
01187     {
01188         if (nZeroArea > 0)
01189         {
01190             Info<< "There are " << nZeroArea
01191                 << " faces with area < " << minArea << '.' << nl << endl;
01192         }
01193         else
01194         {
01195             Info<< "All faces have area > " << minArea << '.' << nl << endl;
01196         }
01197     }
01198 
01199     if (nZeroArea > 0)
01200     {
01201         if (report)
01202         {
01203             WarningIn
01204             (
01205                 "primitiveMeshGeometry::checkFaceArea"
01206                 "(const bool, const scalar, const primitiveMesh&"
01207                 ", const pointField&, const labelList&, labelHashSet*)"
01208             )   << nZeroArea  << " faces with area < " << minArea
01209                 << " found.\n"
01210                 << endl;
01211         }
01212 
01213         return true;
01214     }
01215     else
01216     {
01217         return false;
01218     }
01219 }
01220 
01221 
01222 bool Foam::primitiveMeshGeometry::checkCellDeterminant
01223 (
01224     const bool report,
01225     const scalar warnDet,
01226     const primitiveMesh& mesh,
01227     const vectorField& faceAreas,
01228     const labelList& checkFaces,
01229     const labelList& affectedCells,
01230     labelHashSet* setPtr
01231 )
01232 {
01233     const cellList& cells = mesh.cells();
01234 
01235     scalar minDet = GREAT;
01236     scalar sumDet = 0.0;
01237     label nSumDet = 0;
01238     label nWarnDet = 0;
01239 
01240     forAll(affectedCells, i)
01241     {
01242         const cell& cFaces = cells[affectedCells[i]];
01243 
01244         tensor areaSum(tensor::zero);
01245         scalar magAreaSum = 0;
01246 
01247         forAll(cFaces, cFaceI)
01248         {
01249             label faceI = cFaces[cFaceI];
01250     
01251             scalar magArea = mag(faceAreas[faceI]);
01252 
01253             magAreaSum += magArea;
01254             areaSum += faceAreas[faceI]*(faceAreas[faceI]/magArea);
01255         }
01256 
01257         scalar scaledDet = det(areaSum/magAreaSum)/0.037037037037037;
01258 
01259         minDet = min(minDet, scaledDet);
01260         sumDet += scaledDet;
01261         nSumDet++;
01262 
01263         if (scaledDet < warnDet)
01264         {
01265             if (setPtr)
01266             {
01267                 // Insert all faces of the cell.
01268                 forAll(cFaces, cFaceI)
01269                 {
01270                     label faceI = cFaces[cFaceI];
01271                     setPtr->insert(faceI);
01272                 }
01273             }
01274             nWarnDet++;
01275         }
01276     }
01277     
01278     reduce(minDet, minOp<scalar>());
01279     reduce(sumDet, sumOp<scalar>());
01280     reduce(nSumDet, sumOp<label>());
01281     reduce(nWarnDet, sumOp<label>());
01282 
01283     if (report)
01284     {
01285         if (nSumDet > 0)
01286         {
01287             Info<< "Cell determinant (1 = uniform cube) : average = "
01288                 << sumDet / nSumDet << "  min = " << minDet << endl;
01289         }
01290 
01291         if (nWarnDet > 0)
01292         {
01293             Info<< "There are " << nWarnDet
01294                 << " cells with determinant < " << warnDet << '.' << nl
01295                 << endl;
01296         }
01297         else
01298         {
01299             Info<< "All faces have determinant > " << warnDet << '.' << nl
01300                 << endl;
01301         }
01302     }
01303 
01304     if (nWarnDet > 0)
01305     {
01306         if (report)
01307         {
01308             WarningIn
01309             (
01310                 "primitiveMeshGeometry::checkCellDeterminant"
01311                 "(const bool, const scalar, const primitiveMesh&"
01312                 ", const pointField&, const labelList&, const labelList&"
01313                 ", labelHashSet*)"
01314             )   << nWarnDet << " cells with determinant < " << warnDet
01315                 << " found.\n"
01316                 << endl;
01317         }
01318 
01319         return true;
01320     }
01321     else
01322     {
01323         return false;
01324     }
01325 }
01326 
01327 
01328 bool Foam::primitiveMeshGeometry::checkFaceDotProduct
01329 (
01330     const bool report,
01331     const scalar orthWarn,
01332     const labelList& checkFaces,
01333     labelHashSet* setPtr
01334 ) const
01335 {
01336     return checkFaceDotProduct
01337     (
01338         report,
01339         orthWarn,
01340         mesh_,
01341         cellCentres_,
01342         faceAreas_,
01343         checkFaces,
01344         setPtr
01345     );
01346 }
01347 
01348 
01349 bool Foam::primitiveMeshGeometry::checkFacePyramids
01350 (
01351     const bool report,
01352     const scalar minPyrVol,
01353     const pointField& p,
01354     const labelList& checkFaces,
01355     labelHashSet* setPtr
01356 ) const
01357 {
01358     return checkFacePyramids
01359     (
01360         report,
01361         minPyrVol,
01362         mesh_,
01363         cellCentres_,
01364         p,
01365         checkFaces,
01366         setPtr
01367     );
01368 }
01369 
01370 
01371 bool Foam::primitiveMeshGeometry::checkFaceSkewness
01372 (
01373     const bool report,
01374     const scalar internalSkew,
01375     const scalar boundarySkew,
01376     const labelList& checkFaces,
01377     labelHashSet* setPtr
01378 ) const
01379 {
01380     return checkFaceSkewness
01381     (
01382         report,
01383         internalSkew,
01384         boundarySkew,
01385         mesh_,
01386         cellCentres_,
01387         faceCentres_,
01388         faceAreas_,
01389         checkFaces,
01390         setPtr
01391     );
01392 }
01393 
01394 
01395 bool Foam::primitiveMeshGeometry::checkFaceWeights
01396 (
01397     const bool report,
01398     const scalar warnWeight,
01399     const labelList& checkFaces,
01400     labelHashSet* setPtr
01401 ) const
01402 {
01403     return checkFaceWeights
01404     (
01405         report,
01406         warnWeight,
01407         mesh_,
01408         cellCentres_,
01409         faceCentres_,
01410         faceAreas_,
01411         checkFaces,
01412         setPtr
01413     );
01414 }
01415 
01416 
01417 bool Foam::primitiveMeshGeometry::checkFaceAngles
01418 (
01419     const bool report,
01420     const scalar maxDeg,
01421     const pointField& p,
01422     const labelList& checkFaces,
01423     labelHashSet* setPtr
01424 ) const
01425 {
01426     return checkFaceAngles
01427     (
01428         report,
01429         maxDeg,
01430         mesh_,
01431         faceAreas_,
01432         p,
01433         checkFaces,
01434         setPtr
01435     );
01436 }
01437 
01438 
01439 //bool Foam::primitiveMeshGeometry::checkFaceFlatness
01440 //(
01441 //    const bool report,
01442 //    const scalar warnFlatness,
01443 //    const pointField& p,
01444 //    const labelList& checkFaces,
01445 //    labelHashSet* setPtr
01446 //) const
01447 //{
01448 //    return checkFaceFlatness
01449 //    (
01450 //        report,
01451 //        warnFlatness,
01452 //        mesh_,
01453 //        faceAreas_,
01454 //        faceCentres_,
01455 //        p,
01456 //        checkFaces,
01457 //        setPtr
01458 //    );
01459 //}
01460 
01461 
01462 bool Foam::primitiveMeshGeometry::checkFaceTwist
01463 (
01464     const bool report,
01465     const scalar minTwist,
01466     const pointField& p,
01467     const labelList& checkFaces,
01468     labelHashSet* setPtr
01469 ) const
01470 {
01471     return checkFaceTwist
01472     (
01473         report,
01474         minTwist,
01475         mesh_,
01476         faceAreas_,
01477         faceCentres_,
01478         p,
01479         checkFaces,
01480         setPtr
01481     );
01482 }
01483 
01484 
01485 bool Foam::primitiveMeshGeometry::checkFaceArea
01486 (
01487     const bool report,
01488     const scalar minArea,
01489     const labelList& checkFaces,
01490     labelHashSet* setPtr
01491 ) const
01492 {
01493     return checkFaceArea
01494     (
01495         report,
01496         minArea,
01497         mesh_,
01498         faceAreas_,
01499         checkFaces,
01500         setPtr
01501     );
01502 }
01503 
01504 
01505 bool Foam::primitiveMeshGeometry::checkCellDeterminant
01506 (
01507     const bool report,
01508     const scalar warnDet,
01509     const labelList& checkFaces,
01510     const labelList& affectedCells,
01511     labelHashSet* setPtr
01512 ) const
01513 {
01514     return checkCellDeterminant
01515     (
01516         report,
01517         warnDet,
01518         mesh_,
01519         faceAreas_,
01520         checkFaces,
01521         affectedCells,
01522         setPtr
01523     );
01524 }
01525 
01526 
01527 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines