Fix tutorials: typo in tutorials/viscoelastic/viscoelasticFluidFoam/S-MDCPP/constant...
[OpenFOAM-1.6-ext.git] / src / finiteVolume / interpolation / pointVolInterpolation / pointVolInterpolation.C
blobdac5ab378912b58a6ba7ffbc4958a81d77b5a10a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM 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 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     pointVolInterpolation
28 Description
30 \*---------------------------------------------------------------------------*/
32 #include "pointVolInterpolation.H"
33 #include "fvMesh.H"
34 #include "pointFields.H"
35 #include "volFields.H"
36 #include "emptyFvPatch.H"
37 #include "SubField.H"
38 #include "demandDrivenData.H"
39 #include "globalMeshData.H"
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 namespace Foam
45     defineTypeNameAndDebug(pointVolInterpolation, 0);
49 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
51 void Foam::pointVolInterpolation::makeWeights() const
53     if (volWeightsPtr_)
54     {
55         FatalErrorIn("pointVolInterpolation::makeWeights() const")
56             << "weighting factors already calculated"
57             << abort(FatalError);
58     }
60     if (debug)
61     {
62         Info<< "pointVolInterpolation::makeWeights() : "
63             << "constructing weighting factors"
64             << endl;
65     }
67     const pointField& points = vMesh().points();
68     const labelListList& cellPoints = vMesh().cellPoints();
69     const vectorField& cellCentres = vMesh().cellCentres();
71     // Allocate storage for weighting factors
72     volWeightsPtr_ = new FieldField<Field, scalar>(cellCentres.size());
73     FieldField<Field, scalar>& weightingFactors = *volWeightsPtr_;
75     forAll(weightingFactors, pointi)
76     {
77         weightingFactors.set
78         (
79             pointi,
80             new scalarField(cellPoints[pointi].size())
81         );
82     }
85     // Calculate inverse distances between points and cell centres
86     // and store in weighting factor array
87     forAll(cellCentres, cellI)
88     {
89         const labelList& curCellPoints = cellPoints[cellI];
91         forAll(curCellPoints, cellPointI)
92         {
93             weightingFactors[cellI][cellPointI] = 1.0/
94                 mag
95                 (
96                     cellCentres[cellI] - points[curCellPoints[cellPointI]]
97                 );
98         }
99     }
101     scalarField pointVolSumWeights(cellCentres.size(), 0);
103     // Calculate weighting factors using inverse distance weighting
104     forAll(cellCentres, cellI)
105     {
106         const labelList& curCellPoints = cellPoints[cellI];
108         forAll(curCellPoints, cellPointI)
109         {
110             pointVolSumWeights[cellI] += weightingFactors[cellI][cellPointI];
111         }
112     }
114     forAll(cellCentres, cellI)
115     {
116         const labelList& curCellPoints = cellPoints[cellI];
118         forAll(curCellPoints, cellPointI)
119         {
120             weightingFactors[cellI][cellPointI] /= pointVolSumWeights[cellI];
121         }
122     }
124     if (debug)
125     {
126         Info<< "pointVolInterpolation::makeWeights() : "
127             << "finished constructing weighting factors"
128             << endl;
129     }
133 // Do what is neccessary if the mesh has changed topology
134 void Foam::pointVolInterpolation::clearAddressing() const
136     deleteDemandDrivenData(patchInterpolatorsPtr_);
140 // Do what is neccessary if the mesh has moved
141 void Foam::pointVolInterpolation::clearGeom() const
143     deleteDemandDrivenData(volWeightsPtr_);
147 const Foam::PtrList<Foam::primitivePatchInterpolation>&
148 Foam::pointVolInterpolation::patchInterpolators() const
150     if (!patchInterpolatorsPtr_)
151     {
152         const fvBoundaryMesh& bdry = vMesh().boundary();
154         patchInterpolatorsPtr_ =
155             new PtrList<primitivePatchInterpolation>(bdry.size());
157         forAll (bdry, patchI)
158         {
159             patchInterpolatorsPtr_->set
160             (
161                 patchI,
162                 new primitivePatchInterpolation(bdry[patchI].patch())
163             );
164         }
165     }
167     return *patchInterpolatorsPtr_;
171 // * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * * //
173 Foam::pointVolInterpolation::pointVolInterpolation
175  const pointMesh& pm,
176     const fvMesh& vm
179     pointMesh_(pm),
180     fvMesh_(vm),
181     volWeightsPtr_(NULL),
182     patchInterpolatorsPtr_(NULL)
186 // * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * * //
188 Foam::pointVolInterpolation::~pointVolInterpolation()
190     clearAddressing();
191     clearGeom();
195 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
197 // Return point weights
198 const Foam::FieldField<Foam::Field, Foam::scalar>&
199 Foam::pointVolInterpolation::volWeights() const
201     // If weighting factor array not present then construct
202     if (!volWeightsPtr_)
203     {
204         makeWeights();
205     }
207     return *volWeightsPtr_;
211 // Do what is neccessary if the mesh has moved
212 void Foam::pointVolInterpolation::updateTopology()
214     clearAddressing();
215     clearGeom();
219 // Do what is neccessary if the mesh has moved
220 bool Foam::pointVolInterpolation::movePoints()
222     clearGeom();
224     return true;
228 // ************************************************************************* //