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"
33 #include "volFields.H"
34 #include "pointFields.H"
35 #include "primitiveMesh.H"
36 #include "emptyFvPatch.H"
37 #include "globalMeshData.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 void Foam::pointVolInterpolation::interpolate
44 const GeometricField<Type, pointPatchField, pointMesh>& pf,
45 GeometricField<Type, fvPatchField, volMesh>& vf
50 Info<< "pointVolInterpolation::interpolate("
51 << "const GeometricField<Type, pointPatchField, pointMesh>&, "
52 << "GeometricField<Type, fvPatchField, volMesh>&) : "
53 << "interpolating field from points to cells"
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)
63 vf[cellI] = pTraits<Type>::zero;
65 const labelList& curCellPoints = cellPoints[cellI];
67 forAll(curCellPoints, cellPointI)
70 weights[cellI][cellPointI]*pf[curCellPoints[cellPointI]];
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();
82 // If the patch is empty, skip it
85 bm[patchI].type() != emptyFvPatch::typeName
88 vf.boundaryField()[patchI] =
89 pi[patchI].pointToFaceInterpolate
91 pf.boundaryField()[patchI].patchInternalField()
96 vf.correctBoundaryConditions();
100 Info<< "pointVolInterpolation::interpolate("
101 << "const GeometricField<Type, pointPatchField, pointMesh>&, "
102 << "GeometricField<Type, fvPatchField, volMesh>&) : "
103 << "finished interpolating field from points to cells"
110 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
111 Foam::pointVolInterpolation::interpolate
113 const GeometricField<Type, pointPatchField, pointMesh>& pf
116 // Construct tmp<pointField>
117 tmp<GeometricField<Type, fvPatchField, volMesh> > tvf
119 new GeometricField<Type, fvPatchField, volMesh>
123 "pointVolInterpolate(" + pf.name() + ')',
132 // Perform interpolation
133 interpolate(pf, tvf());
140 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
141 Foam::pointVolInterpolation::interpolate
143 const tmp<GeometricField<Type, pointPatchField, pointMesh> >& tpf
146 // Construct tmp<volField>
147 tmp<GeometricField<Type, fvPatchField, volMesh> > tvf =
154 // ************************************************************************* //