Better bounding on topo change
[foam-extend-3.2.git] / src / dynamicMesh / meshMotion / fvMotionSolver / motionDiffusivity / inverseFaceDistance / inverseFaceDistanceDiffusivity.C
blob9511ecdd8fc4220cc8c1ecdde13472e6dac7bede
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "inverseFaceDistanceDiffusivity.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "HashSet.H"
29 #include "wallPoint.H"
30 #include "MeshWave.H"
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 namespace Foam
36     defineTypeNameAndDebug(inverseFaceDistanceDiffusivity, 0);
38     addToRunTimeSelectionTable
39     (
40         motionDiffusivity,
41         inverseFaceDistanceDiffusivity,
42         Istream
43     );
47 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
49 Foam::inverseFaceDistanceDiffusivity::inverseFaceDistanceDiffusivity
51     const fvMotionSolver& mSolver,
52     Istream& mdData
55     uniformDiffusivity(mSolver, mdData),
56     patchNames_(mdData)
58     correct();
62 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
64 Foam::inverseFaceDistanceDiffusivity::~inverseFaceDistanceDiffusivity()
68 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
70 void Foam::inverseFaceDistanceDiffusivity::correct()
72     const polyMesh& mesh = mSolver().mesh();
73     const polyBoundaryMesh& bdry = mesh.boundaryMesh();
75     labelHashSet patchSet(bdry.size());
77     label nPatchFaces = 0;
79     forAll (patchNames_, i)
80     {
81         label pID = bdry.findPatchID(patchNames_[i]);
83         if (pID > -1)
84         {
85             patchSet.insert(pID);
86             nPatchFaces += bdry[pID].size();
87         }
88     }
90     List<wallPoint> faceDist(nPatchFaces);
91     labelList changedFaces(nPatchFaces);
93     nPatchFaces = 0;
95     forAllConstIter(labelHashSet, patchSet, iter)
96     {
97         const polyPatch& patch = bdry[iter.key()];
99         const vectorField::subField fc = patch.faceCentres();
101         forAll(fc, patchFaceI)
102         {
103             changedFaces[nPatchFaces] = patch.start() + patchFaceI;
105             faceDist[nPatchFaces] = wallPoint(fc[patchFaceI], 0);
107             nPatchFaces++;
108         }
109     }
110     faceDist.setSize(nPatchFaces);
111     changedFaces.setSize(nPatchFaces);
113     MeshWave<wallPoint> waveInfo
114     (
115         mesh,
116         changedFaces,
117         faceDist,
118         mesh.globalData().nTotalCells() // max iterations
119     );
121     const List<wallPoint>& faceInfo = waveInfo.allFaceInfo();
122     const List<wallPoint>& cellInfo = waveInfo.allCellInfo();
124     for (label faceI=0; faceI<mesh.nInternalFaces(); faceI++)
125     {
126         scalar dist = faceInfo[faceI].distSqr();
128         faceDiffusivity_[faceI] = 1.0/sqrt(dist);
129     }
131     forAll(faceDiffusivity_.boundaryField(), patchI)
132     {
133         fvsPatchScalarField& bfld = faceDiffusivity_.boundaryField()[patchI];
135         const unallocLabelList& faceCells = bfld.patch().faceCells();
137         if (patchSet.found(patchI))
138         {
139             forAll(bfld, i)
140             {
141                 scalar dist = cellInfo[faceCells[i]].distSqr();
142                 bfld[i] = 1.0/sqrt(dist);
143             }
144         }
145         else
146         {
147             label start = bfld.patch().patch().start();
149             forAll(bfld, i)
150             {
151                 scalar dist = faceInfo[start+i].distSqr();
152                 bfld[i] = 1.0/sqrt(dist);
153             }
154         }
155     }
159 // ************************************************************************* //