initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / interpolation / volPointInterpolation / volPointInterpolate.C
blob1adc4046c73819ced0e64c00090a69d7286a06f7
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 "volFields.H"
29 #include "pointFields.H"
30 #include "globalPointPatch.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 template<class Type>
40 void volPointInterpolation::interpolateInternalField
42     const GeometricField<Type, fvPatchField, volMesh>& vf,
43     GeometricField<Type, pointPatchField, pointMesh>& pf
44 ) const
46     if (debug)
47     {
48         Info<< "volPointInterpolation::interpolateInternalField("
49             << "const GeometricField<Type, fvPatchField, volMesh>&, "
50             << "GeometricField<Type, pointPatchField, pointMesh>&) : "
51             << "interpolating field from cells to points"
52             << endl;
53     }
55     const labelListList& pointCells = vf.mesh().pointCells();
57     // Multiply volField by weighting factor matrix to create pointField
58     forAll(pointCells, pointi)
59     {
60         const scalarList& pw = pointWeights_[pointi];
61         const labelList& ppc = pointCells[pointi];
63         pf[pointi] = pTraits<Type>::zero;
65         forAll(ppc, pointCelli)
66         {
67             pf[pointi] += pw[pointCelli]*vf[ppc[pointCelli]];
68         }
69     }
73 template<class Type>
74 void volPointInterpolation::interpolate
76     const GeometricField<Type, fvPatchField, volMesh>& vf,
77     GeometricField<Type, pointPatchField, pointMesh>& pf
78 ) const
80     if (debug)
81     {
82         Info<< "volPointInterpolation::interpolate("
83             << "const GeometricField<Type, fvPatchField, volMesh>&, "
84             << "GeometricField<Type, pointPatchField, pointMesh>&) : "
85             << "interpolating field from cells to points"
86             << endl;
87     }
89     interpolateInternalField(vf, pf);
91     // Interpolate to the patches preserving fixed value BCs
92     boundaryInterpolator_.interpolate(vf, pf, false);
96 template<class Type>
97 tmp<GeometricField<Type, pointPatchField, pointMesh> >
98 volPointInterpolation::interpolate
100     const GeometricField<Type, fvPatchField, volMesh>& vf,
101     const wordList& patchFieldTypes
102 ) const
104     wordList types(patchFieldTypes);
106     const pointMesh& pMesh = pointMesh::New(vf.mesh());
108     // If the last patch of the pointBoundaryMesh is the global patch
109     // it must be added to the list of patchField types
110     if
111     (
112         isType<globalPointPatch>
113         (
114             pMesh.boundary()[pMesh.boundary().size() - 1]
115         )
116     )
117     {
118         types.setSize(types.size() + 1);
119         types[types.size()-1] = pMesh.boundary()[types.size()-1].type();
120     }
122     // Construct tmp<pointField>
123     tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf
124     (
125         new GeometricField<Type, pointPatchField, pointMesh>
126         (
127             IOobject
128             (
129                 "volPointInterpolate(" + vf.name() + ')',
130                 vf.instance(),
131                 pMesh.thisDb()
132             ),
133             pMesh,
134             vf.dimensions(),
135             types
136         )
137     );
139     interpolateInternalField(vf, tpf());
141     // Interpolate to the patches overriding fixed value BCs
142     boundaryInterpolator_.interpolate(vf, tpf(), true);
144     return tpf;
148 template<class Type>
149 tmp<GeometricField<Type, pointPatchField, pointMesh> >
150 volPointInterpolation::interpolate
152     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
153     const wordList& patchFieldTypes
154 ) const
156     // Construct tmp<pointField>
157     tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf =
158         interpolate(tvf(), patchFieldTypes);
159     tvf.clear();
160     return tpf;
164 template<class Type>
165 tmp<GeometricField<Type, pointPatchField, pointMesh> >
166 volPointInterpolation::interpolate
168     const GeometricField<Type, fvPatchField, volMesh>& vf
169 ) const
171     const pointMesh& pm = pointMesh::New(vf.mesh());
173     tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf
174     (
175         new GeometricField<Type, pointPatchField, pointMesh>
176         (
177             IOobject
178             (
179                 "volPointInterpolate(" + vf.name() + ')',
180                 vf.instance(),
181                 pm.thisDb()
182             ),
183             pm,
184             vf.dimensions()
185         )
186     );
188     interpolateInternalField(vf, tpf());
189     boundaryInterpolator_.interpolate(vf, tpf(), false);
191     return tpf;
195 template<class Type>
196 tmp<GeometricField<Type, pointPatchField, pointMesh> >
197 volPointInterpolation::interpolate
199     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
200 ) const
202     // Construct tmp<pointField>
203     tmp<GeometricField<Type, pointPatchField, pointMesh> > tpf =
204         interpolate(tvf());
205     tvf.clear();
206     return tpf;
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 } // End namespace Foam
214 // ************************************************************************* //