1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
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
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
30 \*---------------------------------------------------------------------------*/
32 #include "pointVolInterpolation.H"
34 #include "pointFields.H"
35 #include "volFields.H"
36 #include "emptyFvPatch.H"
38 #include "demandDrivenData.H"
39 #include "globalMeshData.H"
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 defineTypeNameAndDebug(pointVolInterpolation, 0);
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 void Foam::pointVolInterpolation::makeWeights() const
55 FatalErrorIn("pointVolInterpolation::makeWeights() const")
56 << "weighting factors already calculated"
62 Info<< "pointVolInterpolation::makeWeights() : "
63 << "constructing weighting factors"
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)
80 new scalarField(cellPoints[pointi].size())
85 // Calculate inverse distances between points and cell centres
86 // and store in weighting factor array
87 forAll(cellCentres, cellI)
89 const labelList& curCellPoints = cellPoints[cellI];
91 forAll(curCellPoints, cellPointI)
93 weightingFactors[cellI][cellPointI] = 1.0/
96 cellCentres[cellI] - points[curCellPoints[cellPointI]]
101 scalarField pointVolSumWeights(cellCentres.size(), 0);
103 // Calculate weighting factors using inverse distance weighting
104 forAll(cellCentres, cellI)
106 const labelList& curCellPoints = cellPoints[cellI];
108 forAll(curCellPoints, cellPointI)
110 pointVolSumWeights[cellI] += weightingFactors[cellI][cellPointI];
114 forAll(cellCentres, cellI)
116 const labelList& curCellPoints = cellPoints[cellI];
118 forAll(curCellPoints, cellPointI)
120 weightingFactors[cellI][cellPointI] /= pointVolSumWeights[cellI];
126 Info<< "pointVolInterpolation::makeWeights() : "
127 << "finished constructing weighting factors"
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_)
152 const fvBoundaryMesh& bdry = vMesh().boundary();
154 patchInterpolatorsPtr_ =
155 new PtrList<primitivePatchInterpolation>(bdry.size());
157 forAll (bdry, patchI)
159 patchInterpolatorsPtr_->set
162 new primitivePatchInterpolation(bdry[patchI].patch())
167 return *patchInterpolatorsPtr_;
171 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
173 Foam::pointVolInterpolation::pointVolInterpolation
181 volWeightsPtr_(NULL),
182 patchInterpolatorsPtr_(NULL)
186 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
188 Foam::pointVolInterpolation::~pointVolInterpolation()
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
207 return *volWeightsPtr_;
211 // Do what is neccessary if the mesh has moved
212 void Foam::pointVolInterpolation::updateTopology()
219 // Do what is neccessary if the mesh has moved
220 bool Foam::pointVolInterpolation::movePoints()
228 // ************************************************************************* //