fix consistancy of gradient on coupled patches
[OpenFOAM-1.6-ext.git] / src / finiteVolume / interpolation / pointVolInterpolation / pointVolInterpolate.C
blob8cac7b99f4ad59b848458c0dcead84674cb69f61
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 "volFields.H"
34 #include "pointFields.H"
35 #include "primitiveMesh.H"
36 #include "emptyFvPatch.H"
37 #include "globalMeshData.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 template<class Type>
42 void Foam::pointVolInterpolation::interpolate
44     const GeometricField<Type, pointPatchField, pointMesh>& pf,
45     GeometricField<Type, fvPatchField, volMesh>& vf
46 ) const
48     if (debug)
49     {
50         Info<< "pointVolInterpolation::interpolate("
51             << "const GeometricField<Type, pointPatchField, pointMesh>&, "
52             << "GeometricField<Type, fvPatchField, volMesh>&) : "
53             << "interpolating field from points to cells"
54             << endl;
55     }
57     const FieldField<Field, scalar>& weights = volWeights();
58     const labelListList& cellPoints = vf.mesh().cellPoints();
60     // Multiply pointField by weighting factor matrix to create volField
61     forAll(cellPoints, cellI)
62     {
63         vf[cellI] = pTraits<Type>::zero;
65         const labelList& curCellPoints = cellPoints[cellI];
67         forAll(curCellPoints, cellPointI)
68         {
69             vf[cellI] +=
70                 weights[cellI][cellPointI]*pf[curCellPoints[cellPointI]];
71         }
72     }
75     // Interpolate patch values: over-ride the internal values for the points
76     // on the patch with the interpolated point values from the faces
77     const fvBoundaryMesh& bm = vMesh().boundary();
79     const PtrList<primitivePatchInterpolation>& pi = patchInterpolators();
80     forAll (bm, patchI)
81     {
82         // If the patch is empty, skip it
83         if
84         (
85             bm[patchI].type() != emptyFvPatch::typeName
86         )
87         {
88             vf.boundaryField()[patchI] =
89                 pi[patchI].pointToFaceInterpolate
90                 (
91                     pf.boundaryField()[patchI].patchInternalField()
92                 );
93         }
94     }
96     vf.correctBoundaryConditions();
98     if (debug)
99     {
100         Info<< "pointVolInterpolation::interpolate("
101             << "const GeometricField<Type, pointPatchField, pointMesh>&, "
102             << "GeometricField<Type, fvPatchField, volMesh>&) : "
103             << "finished interpolating field from points to cells"
104             << endl;
105     }
109 template<class Type>
110 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
111 Foam::pointVolInterpolation::interpolate
113     const GeometricField<Type, pointPatchField, pointMesh>& pf
114 ) const
116     // Construct tmp<pointField>
117     tmp<GeometricField<Type, fvPatchField, volMesh> > tvf
118     (
119         new GeometricField<Type, fvPatchField, volMesh>
120         (
121             IOobject
122             (
123                 "pointVolInterpolate(" + pf.name() + ')',
124                 pf.instance(),
125                 pf.db()
126             ),
127             vMesh(),
128             pf.dimensions()
129         )
130     );
132     // Perform interpolation
133     interpolate(pf, tvf());
135     return tvf;
139 template<class Type>
140 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
141 Foam::pointVolInterpolation::interpolate
143     const tmp<GeometricField<Type, pointPatchField, pointMesh> >& tpf
144 ) const
146     // Construct tmp<volField>
147     tmp<GeometricField<Type, fvPatchField, volMesh> > tvf =
148         interpolate(tpf());
149     tpf.clear();
150     return tvf;
154 // ************************************************************************* //