/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . Application surfaceFeatureExtract Description Extracts and writes surface features to file. All but the basic feature extraction is WIP. \*---------------------------------------------------------------------------*/ #include "surfaceFeatureExtract.H" #include "Time.H" #include "meshTools.H" #include "tensor2D.H" #include "symmTensor2D.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // const Foam::scalar Foam::internalAngleTolerance(80); const Foam::scalar Foam::internalToleranceCosAngle ( cos(degToRad(180 - internalAngleTolerance)) ); const Foam::scalar Foam::externalAngleTolerance(10); const Foam::scalar Foam::externalToleranceCosAngle ( cos(degToRad(180 - externalAngleTolerance)) ); bool Foam::edgesConnected(const edge& e1, const edge& e2) { if ( e1.start() == e2.start() || e1.start() == e2.end() || e1.end() == e2.start() || e1.end() == e2.end() ) { return true; } return false; } Foam::scalar Foam::calcProximityOfFeaturePoints ( const List& hitList, const scalar defaultCellSize ) { scalar minDist = defaultCellSize; for ( label hI1 = 0; hI1 < hitList.size() - 1; ++hI1 ) { const pointIndexHit& pHit1 = hitList[hI1]; if (pHit1.hit()) { for ( label hI2 = hI1 + 1; hI2 < hitList.size(); ++hI2 ) { const pointIndexHit& pHit2 = hitList[hI2]; if (pHit2.hit()) { scalar curDist = mag(pHit1.hitPoint() - pHit2.hitPoint()); minDist = min(curDist, minDist); } } } } return minDist; } Foam::scalar Foam::calcProximityOfFeatureEdges ( const extendedFeatureEdgeMesh& efem, const List& hitList, const scalar defaultCellSize ) { scalar minDist = defaultCellSize; for ( label hI1 = 0; hI1 < hitList.size() - 1; ++hI1 ) { const pointIndexHit& pHit1 = hitList[hI1]; if (pHit1.hit()) { const edge& e1 = efem.edges()[pHit1.index()]; for ( label hI2 = hI1 + 1; hI2 < hitList.size(); ++hI2 ) { const pointIndexHit& pHit2 = hitList[hI2]; if (pHit2.hit()) { const edge& e2 = efem.edges()[pHit2.index()]; // Don't refine if the edges are connected to each other if (!edgesConnected(e1, e2)) { scalar curDist = mag(pHit1.hitPoint() - pHit2.hitPoint()); minDist = min(curDist, minDist); } } } } } return minDist; } void Foam::deleteBox ( const triSurface& surf, const treeBoundBox& bb, const bool removeInside, List& edgeStat ) { forAll(edgeStat, edgeI) { const point eMid = surf.edges()[edgeI].centre(surf.localPoints()); if (removeInside ? bb.contains(eMid) : !bb.contains(eMid)) { edgeStat[edgeI] = surfaceFeatures::NONE; } } } bool Foam::onLine(const point& p, const linePointRef& line) { const point& a = line.start(); const point& b = line.end(); if ( ( p.x() < min(a.x(), b.x()) || p.x() > max(a.x(), b.x()) ) || ( p.y() < min(a.y(), b.y()) || p.y() > max(a.y(), b.y()) ) || ( p.z() < min(a.z(), b.z()) || p.z() > max(a.z(), b.z()) ) ) { return false; } return true; } void Foam::deleteEdges ( const triSurface& surf, const plane& cutPlane, List& edgeStat ) { const pointField& points = surf.points(); const labelList& meshPoints = surf.meshPoints(); forAll(edgeStat, edgeI) { const edge& e = surf.edges()[edgeI]; const point& p0 = points[meshPoints[e.start()]]; const point& p1 = points[meshPoints[e.end()]]; const linePointRef line(p0, p1); // If edge does not intersect the plane, delete. scalar intersect = cutPlane.lineIntersect(line); point featPoint = intersect * (p1 - p0) + p0; if (!onLine(featPoint, line)) { edgeStat[edgeI] = surfaceFeatures::NONE; } } } void Foam::drawHitProblem ( const label fi, const triSurface& surf, const point& start, const point& p, const point& end, const List& hitInfo ) { Info<< nl << "# findLineAll did not hit its own face." << nl << "# fi " << fi << nl << "# start " << start << nl << "# point " << p << nl << "# end " << end << nl << "# hitInfo " << hitInfo << endl; meshTools::writeOBJ(Info, start); meshTools::writeOBJ(Info, p); meshTools::writeOBJ(Info, end); Info<< "l 1 2 3" << endl; meshTools::writeOBJ(Info, surf.points()[surf[fi][0]]); meshTools::writeOBJ(Info, surf.points()[surf[fi][1]]); meshTools::writeOBJ(Info, surf.points()[surf[fi][2]]); Info<< "f 4 5 6" << endl; forAll(hitInfo, hI) { label hFI = hitInfo[hI].index(); meshTools::writeOBJ(Info, surf.points()[surf[hFI][0]]); meshTools::writeOBJ(Info, surf.points()[surf[hFI][1]]); meshTools::writeOBJ(Info, surf.points()[surf[hFI][2]]); Info<< "f " << 3*hI + 7 << " " << 3*hI + 8 << " " << 3*hI + 9 << endl; } } void Foam::unmarkBaffles ( const triSurface& surf, const scalar includedAngle, List& edgeStat ) { scalar minCos = Foam::cos(degToRad(180.0 - includedAngle)); const labelListList& edgeFaces = surf.edgeFaces(); forAll(edgeFaces, edgeI) { const labelList& eFaces = edgeFaces[edgeI]; if (eFaces.size() > 2) { label i0 = eFaces[0]; //const labelledTri& f0 = surf[i0]; const Foam::vector& n0 = surf.faceNormals()[i0]; //Pout<< "edge:" << edgeI << " n0:" << n0 << endl; bool same = true; for (label i = 1; i < eFaces.size(); i++) { //const labelledTri& f = surf[i]; const Foam::vector& n = surf.faceNormals()[eFaces[i]]; //Pout<< " mag(n&n0): " << mag(n&n0) << endl; if (mag(n&n0) < minCos) { same = false; break; } } if (same) { edgeStat[edgeI] = surfaceFeatures::NONE; } } } } Foam::surfaceFeatures::edgeStatus Foam::checkFlatRegionEdge ( const triSurface& surf, const scalar tol, const scalar includedAngle, const label edgeI ) { const edge& e = surf.edges()[edgeI]; const labelList& eFaces = surf.edgeFaces()[edgeI]; // Bin according to normal DynamicList normals(2); DynamicList bins(2); forAll(eFaces, eFacei) { const Foam::vector& n = surf.faceNormals()[eFaces[eFacei]]; // Find the normal in normals label index = -1; forAll(normals, normalI) { if (mag(n&normals[normalI]) > (1-tol)) { index = normalI; break; } } if (index != -1) { bins[index].append(eFacei); } else if (normals.size() >= 2) { // Would be third normal. Mark as feature. //Pout<< "** at edge:" << surf.localPoints()[e[0]] // << surf.localPoints()[e[1]] // << " have normals:" << normals // << " and " << n << endl; return surfaceFeatures::REGION; } else { normals.append(n); bins.append(labelList(1, eFacei)); } } // Check resulting number of bins if (bins.size() == 1) { // Note: should check here whether they are two sets of faces // that are planar or indeed 4 faces al coming together at an edge. //Pout<< "** at edge:" // << surf.localPoints()[e[0]] // << surf.localPoints()[e[1]] // << " have single normal:" << normals[0] // << endl; return surfaceFeatures::NONE; } else { // Two bins. Check if normals make an angle //Pout<< "** at edge:" // << surf.localPoints()[e[0]] // << surf.localPoints()[e[1]] << nl // << " normals:" << normals << nl // << " bins :" << bins << nl // << endl; if (includedAngle >= 0) { scalar minCos = Foam::cos(degToRad(180.0 - includedAngle)); forAll(eFaces, i) { const Foam::vector& ni = surf.faceNormals()[eFaces[i]]; for (label j=i+1; j 0) { regionAndNormal[i] = t.region()+1; } else if (dir == 0) { FatalErrorInFunction << exit(FatalError); } else { regionAndNormal[i] = -(t.region()+1); } } // 2. check against bin1 const labelList& bin1 = bins[1]; labelList regionAndNormal1(bin1.size()); forAll(bin1, i) { const labelledTri& t = surf.localFaces()[eFaces[bin1[i]]]; int dir = t.edgeDirection(e); label myRegionAndNormal; if (dir > 0) { myRegionAndNormal = t.region()+1; } else { myRegionAndNormal = -(t.region()+1); } regionAndNormal1[i] = myRegionAndNormal; label index = findIndex(regionAndNormal, -myRegionAndNormal); if (index == -1) { // Not found. //Pout<< "cannot find region " << myRegionAndNormal // << " in regions " << regionAndNormal << endl; return surfaceFeatures::REGION; } } // Passed all checks, two normal bins with the same contents. //Pout<< "regionAndNormal:" << regionAndNormal << endl; //Pout<< "myRegionAndNormal:" << regionAndNormal1 << endl; return surfaceFeatures::NONE; } } void Foam::writeStats(const extendedFeatureEdgeMesh& fem, Ostream& os) { os << " points : " << fem.points().size() << nl << " of which" << nl << " convex : " << fem.concaveStart() << nl << " concave : " << (fem.mixedStart()-fem.concaveStart()) << nl << " mixed : " << (fem.nonFeatureStart()-fem.mixedStart()) << nl << " non-feature : " << (fem.points().size()-fem.nonFeatureStart()) << nl << " edges : " << fem.edges().size() << nl << " of which" << nl << " external edges : " << fem.internalStart() << nl << " internal edges : " << (fem.flatStart()- fem.internalStart()) << nl << " flat edges : " << (fem.openStart()- fem.flatStart()) << nl << " open edges : " << (fem.multipleStart()- fem.openStart()) << nl << " multiply connected : " << (fem.edges().size()- fem.multipleStart()) << endl; } // ************************************************************************* //