inversePointDistanceDiffusivity.C 5.63 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
OpenFOAM bot's avatar
OpenFOAM bot committed
6 7
     \\/     M anipulation  |
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8 9
    Copyright (C) 2011-2016 OpenFOAM Foundation
    Copyright (C) 2018 OpenCFD Ltd.
10 11 12 13
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

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

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

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

#include "inversePointDistanceDiffusivity.H"
#include "addToRunTimeSelectionTable.H"
31
#include "HashSet.H"
32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53
#include "pointEdgePoint.H"
#include "PointEdgeWave.H"

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

namespace Foam
{
    defineTypeNameAndDebug(inversePointDistanceDiffusivity, 0);

    addToRunTimeSelectionTable
    (
        motionDiffusivity,
        inversePointDistanceDiffusivity,
        Istream
    );
}


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

Foam::inversePointDistanceDiffusivity::inversePointDistanceDiffusivity
(
54
    const fvMesh& mesh,
55 56 57
    Istream& mdData
)
:
58
    uniformDiffusivity(mesh, mdData),
59 60 61 62 63 64 65 66 67 68
    patchNames_(mdData)
{
    correct();
}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

void Foam::inversePointDistanceDiffusivity::correct()
{
69
    const polyBoundaryMesh& bdry = mesh().boundaryMesh();
70

71
    labelHashSet patchSet(bdry.patchSet(patchNames_));
72 73 74

    label nPatchEdges = 0;

75
    for (const label patchi : patchSet)
76
    {
77
        nPatchEdges += bdry[patchi].nEdges();
78 79 80
    }

    // Distance to wall on points and edges.
81 82
    List<pointEdgePoint> pointWallDist(mesh().nPoints());
    List<pointEdgePoint> edgeWallDist(mesh().nEdges());
83

84 85 86
    int dummyTrackData = 0;


87 88 89 90 91 92 93
    {
        // Seeds
        List<pointEdgePoint> seedInfo(nPatchEdges);
        labelList seedPoints(nPatchEdges);

        nPatchEdges = 0;

94
        for (const label patchi : patchSet)
95
        {
96
            const polyPatch& patch = bdry[patchi];
97 98 99

            const labelList& meshPoints = patch.meshPoints();

100
            for (const label pointi : meshPoints)
101
            {
102
                if (!pointWallDist[pointi].valid(dummyTrackData))
103 104 105 106
                {
                    // Not yet seeded
                    seedInfo[nPatchEdges] = pointEdgePoint
                    (
107
                        mesh().points()[pointi],
108 109
                        0.0
                    );
110 111
                    seedPoints[nPatchEdges] = pointi;
                    pointWallDist[pointi] = seedInfo[nPatchEdges];
112 113 114 115 116 117 118 119 120 121 122

                    nPatchEdges++;
                }
            }
        }
        seedInfo.setSize(nPatchEdges);
        seedPoints.setSize(nPatchEdges);

        // Do calculations
        PointEdgeWave<pointEdgePoint> waveInfo
        (
123
            mesh(),
124 125
            seedPoints,
            seedInfo,
126

127 128
            pointWallDist,
            edgeWallDist,
129
            mesh().globalData().nTotalPoints(),// max iterations
130
            dummyTrackData
131 132 133 134
        );
    }


135
    for (label facei=0; facei<mesh().nInternalFaces(); ++facei)
136
    {
137
        const face& f = mesh().faces()[facei];
138 139 140 141 142 143 144 145 146

        scalar dist = 0;

        forAll(f, fp)
        {
            dist += sqrt(pointWallDist[f[fp]].distSqr());
        }
        dist /= f.size();

147
        faceDiffusivity_[facei] = 1.0/dist;
148 149
    }

150

151
    surfaceScalarField::Boundary& faceDiffusivityBf =
152 153
        faceDiffusivity_.boundaryFieldRef();

154
    forAll(faceDiffusivityBf, patchi)
155
    {
156
        fvsPatchScalarField& bfld = faceDiffusivityBf[patchi];
157

158
        if (patchSet.found(patchi))
159
        {
160
            const labelUList& faceCells = bfld.patch().faceCells();
161 162 163

            forAll(bfld, i)
            {
164
                const cell& ownFaces = mesh().cells()[faceCells[i]];
165

166 167 168 169
                labelHashSet cPoints(4*ownFaces.size());

                scalar dist = 0;

170
                forAll(ownFaces, ownFacei)
171
                {
172
                    const face& f = mesh().faces()[ownFaces[ownFacei]];
173 174 175 176 177 178 179 180 181 182

                    forAll(f, fp)
                    {
                        if (cPoints.insert(f[fp]))
                        {
                            dist += sqrt(pointWallDist[f[fp]].distSqr());
                        }
                    }
                }
                dist /= cPoints.size();
183

184 185 186 187 188
                bfld[i] = 1.0/dist;
            }
        }
        else
        {
189
            const label start = bfld.patch().start();
190 191 192

            forAll(bfld, i)
            {
193
                const face& f = mesh().faces()[start+i];
194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210

                scalar dist = 0;

                forAll(f, fp)
                {
                    dist += sqrt(pointWallDist[f[fp]].distSqr());
                }
                dist /= f.size();

                bfld[i] = 1.0/dist;
            }
        }
    }
}


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