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/>.
25 Foam::volPointInterpolation
28 Foam::volPointInterpolation
31 volPointInterpolation.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef volPointInterpolation_H
37 #define volPointInterpolation_H
39 #include "MeshObject.H"
40 #include "pointPatchInterpolation.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 /*---------------------------------------------------------------------------*\
51 Class volPointInterpolation Declaration
52 \*---------------------------------------------------------------------------*/
54 class volPointInterpolation
56 public MeshObject<fvMesh, volPointInterpolation>
60 //- Boundary interpolation engine
61 pointPatchInterpolation boundaryInterpolator_;
63 //- Interpolation scheme weighting factor array
64 // Updated for MeshObject. HJ, 30/Aug/2010
65 mutable scalarListList pointWeights_;
68 // Private member functions
70 //- Construct point weighting factors
71 void makeWeights() const;
73 //- Disallow default bitwise copy construct
74 volPointInterpolation(const volPointInterpolation&);
76 //- Disallow default bitwise assignment
77 void operator=(const volPointInterpolation&);
82 // Declare name of the class and its debug switch
83 TypeName("volPointInterpolation");
88 //- Constructor given fvMesh. pointMesh will be created or
89 // looked up from objectRegistry
90 explicit volPointInterpolation(const fvMesh&);
95 ~volPointInterpolation();
102 const fvMesh& mesh() const
104 return boundaryInterpolator_.mesh();
110 //- Correct weighting factors for moving mesh.
111 // Updated for MeshObject. HJ, 30/Aug/2010
112 virtual bool movePoints() const;
114 //- Update mesh topology using the morph engine
115 // Updated for MeshObject. HJ, 30/Aug/2010
116 virtual bool updateMesh(const mapPolyMesh&) const;
119 // Interpolation functions
121 //- Interpolate internal field from volField to pointField
122 // using inverse distance weighting
124 void interpolateInternalField
126 const GeometricField<Type, fvPatchField, volMesh>&,
127 GeometricField<Type, pointPatchField, pointMesh>&
130 //- Interpolate from volField to pointField
131 // using inverse distance weighting
135 const GeometricField<Type, fvPatchField, volMesh>&,
136 GeometricField<Type, pointPatchField, pointMesh>&
139 //- Interpolate volField using inverse distance weighting
140 // returning pointField with the same patchField types
142 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
144 const GeometricField<Type, fvPatchField, volMesh>&,
145 const wordList& patchFieldTypes
148 //- Interpolate tmp<volField> using inverse distance weighting
149 // returning pointField with the same patchField types
151 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
153 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
154 const wordList& patchFieldTypes
157 //- Interpolate volField using inverse distance weighting
158 // returning pointField
160 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
162 const GeometricField<Type, fvPatchField, volMesh>&
165 //- Interpolate tmp<volField> using inverse distance weighting
166 // returning pointField
168 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
170 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 } // End namespace Foam
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 # include "volPointInterpolate.C"
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 // ************************************************************************* //