1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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 "volPointInterpolation.H"
27 #include "volFields.H"
28 #include "pointFields.H"
29 #include "globalPointPatch.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 void volPointInterpolation::interpolateInternalField
41 const GeometricField<Type, fvPatchField, volMesh>& vf,
42 GeometricField<Type, pointPatchField, pointMesh>& pf
47 Info<< "volPointInterpolation::interpolateInternalField("
48 << "const GeometricField<Type, fvPatchField, volMesh>&, "
49 << "GeometricField<Type, pointPatchField, pointMesh>&) : "
50 << "interpolating field from cells to points"
54 const labelListList& pointCells = vf.mesh().pointCells();
56 // Multiply volField by weighting factor matrix to create pointField
57 forAll(pointCells, pointi)
59 const scalarList& pw = pointWeights_[pointi];
60 const labelList& ppc = pointCells[pointi];
62 pf[pointi] = pTraits<Type>::zero;
64 forAll(ppc, pointCelli)
66 pf[pointi] += pw[pointCelli]*vf[ppc[pointCelli]];
73 void volPointInterpolation::interpolate
75 const GeometricField<Type, fvPatchField, volMesh>& vf,
76 GeometricField<Type, pointPatchField, pointMesh>& pf
81 Info<< "volPointInterpolation::interpolate("
82 << "const GeometricField<Type, fvPatchField, volMesh>&, "
83 << "GeometricField<Type, pointPatchField, pointMesh>&) : "
84 << "interpolating field from cells to points"
88 interpolateInternalField(vf, pf);
90 // Interpolate to the patches preserving fixed value BCs
91 boundaryInterpolator_.interpolate(vf, pf, false);
96 tmp<GeometricField<Type, pointPatchField, pointMesh> >
97 volPointInterpolation::interpolate
99 const GeometricField<Type, fvPatchField, volMesh>& vf,
100 const wordList& patchFieldTypes
103 wordList types(patchFieldTypes);
105 const pointMesh& pMesh = pointMesh::New(vf.mesh());
107 // If the last patch of the pointBoundaryMesh is the global patch
108 // it must be added to the list of patchField types
111 isType<globalPointPatch>
113 pMesh.boundary()[pMesh.boundary().size() - 1]
117 types.setSize(types.size() + 1);
118 types[types.size()-1] = pMesh.boundary()[types.size()-1].type();
121 // Construct tmp<pointField>
122 tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf
124 new GeometricField<Type, pointPatchField, pointMesh>
128 "volPointInterpolate(" + vf.name() + ')',
138 interpolateInternalField(vf, tpf());
140 // Interpolate to the patches overriding fixed value BCs
141 boundaryInterpolator_.interpolate(vf, tpf(), true);
148 tmp<GeometricField<Type, pointPatchField, pointMesh> >
149 volPointInterpolation::interpolate
151 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
152 const wordList& patchFieldTypes
155 // Construct tmp<pointField>
156 tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf =
157 interpolate(tvf(), patchFieldTypes);
164 tmp<GeometricField<Type, pointPatchField, pointMesh> >
165 volPointInterpolation::interpolate
167 const GeometricField<Type, fvPatchField, volMesh>& vf
170 const pointMesh& pm = pointMesh::New(vf.mesh());
172 tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf
174 new GeometricField<Type, pointPatchField, pointMesh>
178 "volPointInterpolate(" + vf.name() + ')',
187 interpolateInternalField(vf, tpf());
188 boundaryInterpolator_.interpolate(vf, tpf(), false);
195 tmp<GeometricField<Type, pointPatchField, pointMesh> >
196 volPointInterpolation::interpolate
198 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
201 // Construct tmp<pointField>
202 tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf =
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 } // End namespace Foam
213 // ************************************************************************* //