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

geomCellLooper.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 "geomCellLooper.H"
00027 #include <OpenFOAM/polyMesh.H>
00028 #include <OpenFOAM/DynamicList.H>
00029 #include <OpenFOAM/plane.H>
00030 #include <meshTools/meshTools.H>
00031 #include <OpenFOAM/SortableList.H>
00032 #include <meshTools/triSurfaceTools.H>
00033 #include <OpenFOAM/HashSet.H>
00034 #include <OpenFOAM/ListOps.H>
00035 #include <OpenFOAM/transform.H>
00036 
00037 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00038 
00039 
00040 // Extension factor of edges to make sure we catch intersections through
00041 // edge endpoints
00042 const Foam::scalar Foam::geomCellLooper::pointEqualTol_ = 1E-3;
00043 
00044 
00045 // Snap cuts through edges onto edge endpoints. Fraction of edge length.
00046 Foam::scalar Foam::geomCellLooper::snapTol_ = 0.1;
00047 
00048 
00049 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00050 
00051 namespace Foam
00052 {
00053 
00054 defineTypeNameAndDebug(geomCellLooper, 0);
00055 addToRunTimeSelectionTable(cellLooper, geomCellLooper, word);
00056 
00057 }
00058 
00059 
00060 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00061 
00062 Foam::scalar Foam::geomCellLooper::minEdgeLen(const label vertI) const
00063 {
00064     scalar minLen = GREAT;
00065 
00066     const labelList& pEdges = mesh().pointEdges()[vertI];
00067 
00068     forAll(pEdges, pEdgeI)
00069     {
00070         const edge& e = mesh().edges()[pEdges[pEdgeI]];
00071 
00072         minLen = min(minLen, e.mag(mesh().points()));
00073     }
00074     return minLen;
00075 }
00076 
00077 
00078 // Cut edge with plane. Return true and set weight to fraction between
00079 // edge-start and edge-end
00080 bool Foam::geomCellLooper::cutEdge
00081 (
00082     const plane& cutPlane,
00083     const label edgeI,
00084     scalar& weight
00085 ) const
00086 {
00087     const pointField& pts = mesh().points();
00088 
00089     const edge& e = mesh().edges()[edgeI];
00090 
00091     scalar s = cutPlane.normalIntersect(pts[e.start()], e.vec(pts));
00092 
00093     if ((s > -pointEqualTol_) && (s < 1 + pointEqualTol_))
00094     {
00095         weight = s;
00096 
00097         return true;
00098     }
00099     else
00100     {
00101         // Make sure we don't use this value
00102         weight = -GREAT;
00103 
00104         return false;
00105     }
00106 }
00107 
00108 
00109 // If edge close enough to endpoint snap to endpoint.
00110 Foam::label Foam::geomCellLooper::snapToVert
00111 (
00112     const scalar tol,
00113     const label edgeI,
00114     const scalar weight
00115 ) const
00116 {
00117     const edge& e = mesh().edges()[edgeI];
00118 
00119     if (weight < tol)
00120     {
00121         return e.start();
00122     }
00123     else if (weight > (1-tol))
00124     {
00125         return e.end();
00126     }
00127     else
00128     {
00129         return -1;
00130     }
00131 }
00132 
00133 
00134 void Foam::geomCellLooper::getBase(const vector& n, vector& e0, vector& e1)
00135  const
00136 {
00137     // Guess for vector normal to n.
00138     vector base(1,0,0);
00139 
00140     scalar nComp = n & base;
00141 
00142     if (mag(nComp) > 0.8)
00143     {
00144         // Was bad guess. Try with different vector.
00145 
00146         base.x() = 0;
00147         base.y() = 1;
00148 
00149         nComp = n & base;
00150 
00151         if (mag(nComp) > 0.8)
00152         {
00153             base.y() = 0;
00154             base.z() = 1;
00155 
00156             nComp = n & base;
00157         }
00158     }
00159 
00160 
00161     // Use component normal to n as base vector.
00162     e0 = base - nComp*n;
00163 
00164     e0 /= mag(e0) + VSMALL;
00165 
00166     e1 = n ^ e0;
00167 
00168     //Pout<< "Coord system:" << endl
00169     //    << "    n  : " << n << ' ' << mag(n) << endl
00170     //    << "    e0 : " << e0 << ' ' << mag(e0) << endl
00171     //    << "    e1 : " << e1 << ' ' << mag(e1) << endl
00172     //    << endl;
00173 }
00174 
00175 
00176 // Return true if the cut edge at loop[index] is preceded by cuts through
00177 // the edge end points.
00178 bool Foam::geomCellLooper::edgeEndsCut
00179 (
00180     const labelList& loop,
00181     const label index
00182 ) const
00183 {
00184     label edgeI = getEdge(loop[index]);
00185 
00186     const edge& e = mesh().edges()[edgeI];
00187 
00188     const label prevCut = loop[loop.rcIndex(index)];
00189     const label nextCut = loop[loop.fcIndex(index)];
00190 
00191     if (!isEdge(prevCut) && !isEdge(nextCut))
00192     {
00193         // Cuts before and after are both vertices. Check if both
00194         // the edge endpoints
00195         label v0 = getVertex(prevCut);
00196         label v1 = getVertex(nextCut);
00197 
00198         if
00199         (
00200             (v0 == e[0] && v1 == e[1])
00201          || (v0 == e[1] && v1 == e[0])
00202         )
00203         {
00204             return true;
00205         }
00206     }
00207     return false;
00208 }
00209 
00210 
00211 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00212 
00213 // Construct from components
00214 Foam::geomCellLooper::geomCellLooper(const polyMesh& mesh)
00215 :
00216     cellLooper(mesh)
00217 {}
00218 
00219 
00220 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00221 
00222 Foam::geomCellLooper::~geomCellLooper()
00223 {}
00224 
00225 
00226 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00227 
00228 bool Foam::geomCellLooper::cut
00229 (
00230     const vector& refDir,
00231     const label cellI,
00232     const boolList& vertIsCut,
00233     const boolList& edgeIsCut,
00234     const scalarField& edgeWeight,
00235 
00236     labelList& loop,
00237     scalarField& loopWeights
00238 ) const
00239 {
00240     // Cut through cell centre normal to refDir.
00241     return cut
00242     (
00243         plane(mesh().cellCentres()[cellI], refDir),
00244         cellI,
00245         vertIsCut,
00246         edgeIsCut,
00247         edgeWeight,
00248         loop,
00249         loopWeights
00250     );
00251 }
00252 
00253 
00254 bool Foam::geomCellLooper::cut
00255 (
00256     const plane& cutPlane,
00257     const label cellI,
00258     const boolList&,
00259     const boolList&,
00260     const scalarField&,
00261 
00262     labelList& loop,
00263     scalarField& loopWeights
00264 ) const
00265 {
00266     const pointField& points = mesh().points();
00267     const edgeList& edges = mesh().edges();
00268 
00269     // Find all cuts through edges.
00270     // Special cases:
00271     // - edge cut close to endpoint. Snap to endpoint.
00272     // - edge endpoint close to plane (but edge not cut). Snap to endpoint.
00273     // - both endpoints close to plane. Edge parallel to plane. No need to
00274     //   cut to edge.
00275     // Note: any snap-to-point check must be based on all edges using a point
00276     // since otherwise cut through close to point but on neighbouring edge
00277     // might not be snapped.
00278 
00279     // Size overly big.
00280     label nEstCuts = 2*mesh().cells()[cellI].size();
00281 
00282     DynamicList<label> localLoop(nEstCuts);
00283     DynamicList<scalar> localLoopWeights(nEstCuts);
00284 
00285     // Points checked. Used to make sure we don't cut edge and edge endpoints
00286     // at the same time.
00287     labelHashSet checkedPoints(nEstCuts);
00288 
00289     const labelList& cellEdges = mesh().cellEdges()[cellI];
00290 
00291     forAll(cellEdges, i)
00292     {
00293         label edgeI = cellEdges[i];
00294 
00295         const edge& e = edges[edgeI];
00296 
00297         bool useStart = false;
00298 
00299         bool useEnd = false;
00300 
00301         //
00302         // Check distance of endpoints to cutPlane
00303         //
00304 
00305         if (!checkedPoints.found(e.start()))
00306         {
00307             checkedPoints.insert(e.start());
00308 
00309             scalar typStartLen = pointEqualTol_ * minEdgeLen(e.start());
00310 
00311             // Check distance of startPt to plane.
00312             if (cutPlane.distance(points[e.start()]) < typStartLen)
00313             {
00314                 // Use point.
00315                 localLoop.append(vertToEVert(e.start()));
00316                 localLoopWeights.append(-GREAT);
00317 
00318                 useStart = true;
00319             }
00320         }
00321         if (!checkedPoints.found(e.end()))
00322         {
00323             checkedPoints.insert(e.end());
00324 
00325             scalar typEndLen = pointEqualTol_ * minEdgeLen(e.end());
00326 
00327             // Check distance of endPt to plane.
00328             if (cutPlane.distance(points[e.end()]) < typEndLen)
00329             {
00330                 // Use point.
00331                 localLoop.append(vertToEVert(e.end()));
00332                 localLoopWeights.append(-GREAT);
00333 
00334                 useEnd = true;
00335             }
00336         }
00337 
00338         //
00339         // Check cut of edge
00340         //
00341 
00342         if (!useEnd && !useStart)
00343         {
00344             // Edge end points not close to plane so edge not close to
00345             // plane. Cut edge.
00346             scalar cutWeight;
00347 
00348             if (cutEdge(cutPlane, edgeI, cutWeight))
00349             {
00350                 // Snap edge cut to vertex.
00351                 label cutVertI = snapToVert(snapTol_, edgeI, cutWeight);
00352 
00353                 if (cutVertI == -1)
00354                 {
00355                     // Proper cut of edge
00356                     localLoop.append(edgeToEVert(edgeI));
00357                     localLoopWeights.append(cutWeight);
00358                 }
00359                 else
00360                 {
00361                     // Cut through edge end point. Might be duplicate
00362                     // since connected edge might also be snapped to same
00363                     // endpoint so only insert if unique.
00364                     label cut = vertToEVert(cutVertI);
00365 
00366                     if (findIndex(localLoop, cut) == -1)
00367                     {
00368                         localLoop.append(vertToEVert(cutVertI));
00369                         localLoopWeights.append(-GREAT);
00370                     }
00371                 }
00372             }
00373         }
00374     }
00375 
00376     if (localLoop.size() <= 2)
00377     {
00378         return false;
00379     }
00380 
00381     localLoop.shrink();
00382     localLoopWeights.shrink();
00383 
00384 
00385     // Get points on loop and centre of loop
00386     pointField loopPoints(localLoop.size());
00387     point ctr(vector::zero);
00388 
00389     forAll(localLoop, i)
00390     {
00391         loopPoints[i] = coord(localLoop[i], localLoopWeights[i]);
00392         ctr += loopPoints[i];
00393     }
00394     ctr /= localLoop.size();
00395 
00396 
00397     // Get base vectors of coordinate system normal to refDir
00398     vector e0, e1;
00399     getBase(cutPlane.normal(), e0, e1);
00400 
00401 
00402     // Get sorted angles from point on loop to centre of loop.
00403     SortableList<scalar> sortedAngles(localLoop.size());
00404 
00405     forAll(sortedAngles, i)
00406     {
00407         vector toCtr(loopPoints[i] - ctr);
00408         toCtr /= mag(toCtr);
00409 
00410         sortedAngles[i] = pseudoAngle(e0, e1, toCtr);
00411     }
00412     sortedAngles.sort();
00413 
00414     loop.setSize(localLoop.size());
00415     loopWeights.setSize(loop.size());
00416 
00417 
00418     // Fill loop and loopweights according to sorted angles
00419 
00420     const labelList& indices = sortedAngles.indices();
00421 
00422     forAll(indices, i)
00423     {
00424         loop[i] = localLoop[indices[i]];
00425         loopWeights[i] = localLoopWeights[indices[i]];
00426     }
00427 
00428 
00429     // Check for cut edges along already cut vertices.
00430     bool filterLoop = false;
00431 
00432     forAll(loop, i)
00433     {
00434         label cut = loop[i];
00435 
00436         if (isEdge(cut) && edgeEndsCut(loop, i))
00437         {
00438             filterLoop = true;
00439         }
00440     }
00441 
00442     if (filterLoop)
00443     {
00444         // Found edges in loop where both end points are also in loop. This
00445         // is superfluous information so we can remove it.
00446 
00447         labelList filteredLoop(loop.size());
00448         scalarField filteredLoopWeights(loopWeights.size());
00449         label filterI = 0;
00450 
00451         forAll(loop, i)
00452         {
00453             label cut = loop[i];
00454 
00455             if (isEdge(cut) && edgeEndsCut(loop, i))
00456             {
00457                 // Superfluous cut. Edge end points cut and edge itself as well.
00458             }
00459             else
00460             {
00461                 filteredLoop[filterI] = loop[i];
00462                 filteredLoopWeights[filterI] = loopWeights[i];
00463                 filterI++;
00464             }
00465         }
00466         filteredLoop.setSize(filterI);
00467         filteredLoopWeights.setSize(filterI);
00468 
00469         loop.transfer(filteredLoop);
00470         loopWeights.transfer(filteredLoopWeights);
00471     }
00472 
00473     if (debug&2)
00474     {
00475         Pout<< "cell:" << cellI << endl;
00476         forAll(loop, i)
00477         {
00478             Pout<< "At angle:" << sortedAngles[i] << endl
00479                 << "    cut:";
00480 
00481             writeCut(Pout, loop[i], loopWeights[i]);
00482 
00483             Pout<< "  coord:" << coord(loop[i], loopWeights[i]);
00484 
00485             Pout<< endl;
00486         }
00487     }
00488 
00489     return true;
00490 }
00491 
00492 
00493 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines