cellPointWeight.C 6.85 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
6
     \\/     M anipulation  |
OpenFOAM bot's avatar
OpenFOAM bot committed
7
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8
    Copyright (C) 2011-2017 OpenFOAM Foundation
9 10 11 12
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

13 14 15 16
    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.
17 18 19 20 21 22 23

    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
24
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
25 26 27 28 29

\*---------------------------------------------------------------------------*/

#include "cellPointWeight.H"
#include "polyMesh.H"
30
#include "polyMeshTetDecomposition.H"
31 32 33 34

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

int Foam::cellPointWeight::debug(debug::debugSwitch("cellPointWeight", 0));
35

36 37
Foam::scalar Foam::cellPointWeight::tol(SMALL);

Henry Weller's avatar
Henry Weller committed
38
// * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
39 40 41 42 43

void Foam::cellPointWeight::findTetrahedron
(
    const polyMesh& mesh,
    const vector& position,
44
    const label celli
45 46 47 48
)
{
    if (debug)
    {
49
        Pout<< nl << "Foam::cellPointWeight::findTetrahedron" << nl
50
            << "position = " << position << nl
51
            << "celli = " << celli << endl;
52 53
    }

54 55 56
    List<tetIndices> cellTets = polyMeshTetDecomposition::cellTetIndices
    (
        mesh,
57
        celli
58
    );
59

60
    const scalar cellVolume = mesh.cellVolumes()[celli];
61

62
    forAll(cellTets, tetI)
63
    {
64 65 66
        const tetIndices& tetIs = cellTets[tetI];

        // Barycentric coordinates of the position
67
        scalar det = tetIs.tet(mesh).pointToBarycentric(position, weights_);
68

69
        if (mag(det/cellVolume) > tol)
70
        {
71 72 73 74 75 76 77 78 79 80 81 82
            const scalar& u = weights_[0];
            const scalar& v = weights_[1];
            const scalar& w = weights_[2];

            if
            (
                (u + tol > 0)
             && (v + tol > 0)
             && (w + tol > 0)
             && (u + v + w < 1 + tol)
            )
            {
83 84

                faceVertices_ = tetIs.faceTriIs(mesh);
85

86 87 88 89
                return;
            }
        }
    }
90

91 92
    // A suitable point in a tetrahedron was not found, find the
    // nearest.
93

94
    scalar minNearDist = VGREAT;
95

96
    label nearestTetI = -1;
97

98 99 100 101 102
    forAll(cellTets, tetI)
    {
        const tetIndices& tetIs = cellTets[tetI];

        scalar nearDist = tetIs.tet(mesh).nearestPoint(position).distance();
103

104 105 106 107 108
        if (nearDist < minNearDist)
        {
            minNearDist = nearDist;

            nearestTetI = tetI;
109 110 111
        }
    }

112 113 114
    if (debug)
    {
        Pout<< "cellPointWeight::findTetrahedron" << nl
115 116 117
            << "    Tetrahedron search failed; using closest tet to point "
            << position << nl
            << "    cell: "
118
            << celli << nl
119
            << endl;
120
    }
121 122


123 124 125 126 127
    const tetIndices& tetIs = cellTets[nearestTetI];

    // Barycentric coordinates of the position, ignoring if the
    // determinant is suitable.  If not, the return from barycentric
    // to weights_ is safe.
128
    weights_ = tetIs.tet(mesh).pointToBarycentric(position);
129

130
    faceVertices_ = tetIs.faceTriIs(mesh);
131 132 133 134 135 136 137
}


void Foam::cellPointWeight::findTriangle
(
    const polyMesh& mesh,
    const vector& position,
138
    const label facei
139 140 141 142 143 144
)
{
    if (debug)
    {
        Pout<< "\nbool Foam::cellPointWeight::findTriangle" << nl
            << "position = " << position << nl
145
            << "facei = " << facei << endl;
146 147
    }

148 149 150
    List<tetIndices> faceTets = polyMeshTetDecomposition::faceTetIndices
    (
        mesh,
151 152
        facei,
        mesh.faceOwner()[facei]
153
    );
154

155
    const scalar faceAreaSqr = magSqr(mesh.faceAreas()[facei]);
156

157
    forAll(faceTets, tetI)
158
    {
159 160 161
        const tetIndices& tetIs = faceTets[tetI];

        // Barycentric coordinates of the position
162 163 164
        barycentric2D triWeights;
        const scalar det =
            tetIs.faceTri(mesh).pointToBarycentric(position, triWeights);
165 166 167 168 169 170 171 172 173 174 175 176

        if (0.25*mag(det)/faceAreaSqr > tol)
        {
            const scalar& u = triWeights[0];
            const scalar& v = triWeights[1];

            if
            (
                (u + tol > 0)
             && (v + tol > 0)
             && (u + v < 1 + tol)
            )
177
            {
178 179 180 181 182
                // Weight[0] is for the cell centre.
                weights_[0] = 0;
                weights_[1] = triWeights[0];
                weights_[2] = triWeights[1];
                weights_[3] = triWeights[2];
183

184
                faceVertices_ = tetIs.faceTriIs(mesh);
185 186 187 188

                return;
            }
        }
189 190 191 192 193 194 195 196 197 198 199 200 201
    }

    // A suitable point in a triangle was not found, find the nearest.

    scalar minNearDist = VGREAT;

    label nearestTetI = -1;

    forAll(faceTets, tetI)
    {
        const tetIndices& tetIs = faceTets[tetI];

        scalar nearDist = tetIs.faceTri(mesh).nearestPoint(position).distance();
202

203 204 205 206 207 208
        if (nearDist < minNearDist)
        {
            minNearDist = nearDist;

            nearestTetI = tetI;
        }
209 210
    }

211 212
    if (debug)
    {
213 214 215 216
        Pout<< "cellPointWeight::findTriangle" << nl
            << "    Triangle search failed; using closest tri to point "
            << position << nl
            << "    face: "
217
            << facei << nl
218
            << endl;
219
    }
220

221 222 223 224 225 226
    const tetIndices& tetIs = faceTets[nearestTetI];

    // Barycentric coordinates of the position, ignoring if the
    // determinant is suitable.  If not, the return from barycentric
    // to triWeights is safe.

227 228
    const barycentric2D triWeights =
        tetIs.faceTri(mesh).pointToBarycentric(position);
229 230 231 232 233 234

    // Weight[0] is for the cell centre.
    weights_[0] = 0;
    weights_[1] = triWeights[0];
    weights_[2] = triWeights[1];
    weights_[3] = triWeights[2];
235

236
    faceVertices_ = tetIs.faceTriIs(mesh);
237 238 239 240 241 242 243 244 245
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::cellPointWeight::cellPointWeight
(
    const polyMesh& mesh,
    const vector& position,
246 247
    const label celli,
    const label facei
248 249
)
:
250
    celli_(celli)
251
{
252
    if (facei < 0)
253 254
    {
        // Face data not supplied
255
        findTetrahedron(mesh, position, celli);
256 257 258 259
    }
    else
    {
        // Face data supplied
260
        findTriangle(mesh, position, facei);
261 262 263 264 265
    }
}


// ************************************************************************* //