HOTFIX: Adjust Windows build instructions repo
[foam-extend-3.2.git] / src / finiteVolume / interpolation / pointVolInterpolation / pointVolInterpolate.C
blobf76864652e3c4aa34be08dc61e9db64e7a0ead87
1 /*---------------------------------------------------------------------------*\
2   =========                 |
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 -------------------------------------------------------------------------------
8 License
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 Class
25     pointVolInterpolation
27 Description
29 \*---------------------------------------------------------------------------*/
31 #include "pointVolInterpolation.H"
32 #include "volFields.H"
33 #include "pointFields.H"
34 #include "primitiveMesh.H"
35 #include "emptyFvPatch.H"
36 #include "globalMeshData.H"
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 template<class Type>
41 void Foam::pointVolInterpolation::interpolate
43     const GeometricField<Type, pointPatchField, pointMesh>& pf,
44     GeometricField<Type, fvPatchField, volMesh>& vf
45 ) const
47     if (debug)
48     {
49         Info<< "pointVolInterpolation::interpolate("
50             << "const GeometricField<Type, pointPatchField, pointMesh>&, "
51             << "GeometricField<Type, fvPatchField, volMesh>&) : "
52             << "interpolating field from points to cells"
53             << endl;
54     }
56     const FieldField<Field, scalar>& weights = volWeights();
57     const labelListList& cellPoints = vf.mesh().cellPoints();
59     // Multiply pointField by weighting factor matrix to create volField
60     forAll(cellPoints, cellI)
61     {
62         vf[cellI] = pTraits<Type>::zero;
64         const labelList& curCellPoints = cellPoints[cellI];
66         forAll(curCellPoints, cellPointI)
67         {
68             vf[cellI] +=
69                 weights[cellI][cellPointI]*pf[curCellPoints[cellPointI]];
70         }
71     }
74     // Interpolate patch values: over-ride the internal values for the points
75     // on the patch with the interpolated point values from the faces
76     const fvBoundaryMesh& bm = vMesh().boundary();
78     const PtrList<primitivePatchInterpolation>& pi = patchInterpolators();
79     forAll (bm, patchI)
80     {
81         // If the patch is empty, skip it
82         if
83         (
84             bm[patchI].type() != emptyFvPatch::typeName
85         )
86         {
87             vf.boundaryField()[patchI] =
88                 pi[patchI].pointToFaceInterpolate
89                 (
90                     pf.boundaryField()[patchI].patchInternalField()
91                 );
92         }
93     }
95     vf.correctBoundaryConditions();
97     if (debug)
98     {
99         Info<< "pointVolInterpolation::interpolate("
100             << "const GeometricField<Type, pointPatchField, pointMesh>&, "
101             << "GeometricField<Type, fvPatchField, volMesh>&) : "
102             << "finished interpolating field from points to cells"
103             << endl;
104     }
108 template<class Type>
109 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
110 Foam::pointVolInterpolation::interpolate
112     const GeometricField<Type, pointPatchField, pointMesh>& pf
113 ) const
115     // Construct tmp<pointField>
116     tmp<GeometricField<Type, fvPatchField, volMesh> > tvf
117     (
118         new GeometricField<Type, fvPatchField, volMesh>
119         (
120             IOobject
121             (
122                 "pointVolInterpolate(" + pf.name() + ')',
123                 pf.instance(),
124                 pf.db()
125             ),
126             vMesh(),
127             pf.dimensions()
128         )
129     );
131     // Perform interpolation
132     interpolate(pf, tvf());
134     return tvf;
138 template<class Type>
139 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
140 Foam::pointVolInterpolation::interpolate
142     const tmp<GeometricField<Type, pointPatchField, pointMesh> >& tpf
143 ) const
145     // Construct tmp<volField>
146     tmp<GeometricField<Type, fvPatchField, volMesh> > tvf =
147         interpolate(tpf());
148     tpf.clear();
149     return tvf;
153 // ************************************************************************* //