initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / interpolation / volPointInterpolation / volPointInterpolation.C
blobf94ff4456ab3dc085da3c107def0d455411e21ec
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 \*---------------------------------------------------------------------------*/
27 #include "volPointInterpolation.H"
28 #include "fvMesh.H"
29 #include "volFields.H"
30 #include "pointFields.H"
31 #include "demandDrivenData.H"
32 #include "coupledPointPatchFields.H"
33 #include "pointConstraint.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace Foam
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 defineTypeNameAndDebug(volPointInterpolation, 0);
45 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
47 void volPointInterpolation::makeWeights()
49     if (debug)
50     {
51         Info<< "volPointInterpolation::makeWeights() : "
52             << "constructing weighting factors"
53             << endl;
54     }
56     const pointField& points = mesh().points();
57     const labelListList& pointCells = mesh().pointCells();
58     const vectorField& cellCentres = mesh().cellCentres();
60     // Allocate storage for weighting factors
61     pointWeights_.clear();
62     pointWeights_.setSize(points.size());
64     forAll(pointWeights_, pointi)
65     {
66         pointWeights_[pointi].setSize(pointCells[pointi].size());
67     }
69     pointScalarField sumWeights
70     (
71         IOobject
72         (
73             "volPointSumWeights",
74             mesh().polyMesh::instance(),
75             mesh()
76         ),
77         pointMesh::New(mesh()),
78         dimensionedScalar("zero", dimless, 0)
79     );
81     // Calculate inverse distances between cell centres and points
82     // and store in weighting factor array
83     forAll(points, pointi)
84     {
85         scalarList& pw = pointWeights_[pointi];
86         const labelList& pcp = pointCells[pointi];
88         forAll(pcp, pointCelli)
89         {
90             pw[pointCelli] = 
91                 1.0/mag(points[pointi] - cellCentres[pcp[pointCelli]]);
93             sumWeights[pointi] += pw[pointCelli];
94         }
95     }
97     forAll(sumWeights.boundaryField(), patchi)
98     {
99         if (sumWeights.boundaryField()[patchi].coupled())
100         {
101             refCast<coupledPointPatchScalarField>
102                 (sumWeights.boundaryField()[patchi]).initSwapAdd
103                 (
104                     sumWeights.internalField()
105                 );
106         }
107     }
109     forAll(sumWeights.boundaryField(), patchi)
110     {
111         if (sumWeights.boundaryField()[patchi].coupled())
112         {
113             refCast<coupledPointPatchScalarField>
114             (sumWeights.boundaryField()[patchi]).swapAdd
115             (
116                 sumWeights.internalField()
117             );
118         }
119     }
121     forAll(points, pointi)
122     {
123         scalarList& pw = pointWeights_[pointi];
125         forAll(pw, pointCelli)
126         {
127             pw[pointCelli] /= sumWeights[pointi];
128         }
129     }
131     if (debug)
132     {
133         Info<< "volPointInterpolation::makeWeights() : "
134             << "finished constructing weighting factors"
135             << endl;
136     }
140 // * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * * //
142 volPointInterpolation::volPointInterpolation(const fvMesh& vm)
144     MeshObject<fvMesh, volPointInterpolation>(vm),
145     boundaryInterpolator_(vm)
147     updateMesh();
151 // * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * * //
153 volPointInterpolation::~volPointInterpolation()
157 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
159 void volPointInterpolation::updateMesh()
161     makeWeights();
162     boundaryInterpolator_.updateMesh();
166 bool volPointInterpolation::movePoints()
168     makeWeights();
169     boundaryInterpolator_.movePoints();
171     return true;
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 } // End namespace Foam
179 // ************************************************************************* //