1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, 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 scalarListList pointWeights_;
67 // Private member functions
69 //- Construct point weighting factors
72 //- Disallow default bitwise copy construct
73 volPointInterpolation(const volPointInterpolation&);
75 //- Disallow default bitwise assignment
76 void operator=(const volPointInterpolation&);
81 // Declare name of the class and its debug switch
82 ClassName("volPointInterpolation");
87 //- Constructor given fvMesh and pointMesh.
88 explicit volPointInterpolation(const fvMesh&);
93 ~volPointInterpolation();
100 const fvMesh& mesh() const
102 return boundaryInterpolator_.mesh();
108 //- Update mesh topology using the morph engine
111 //- Correct weighting factors for moving mesh.
115 // Interpolation functions
117 //- Interpolate internal field from volField to pointField
118 // using inverse distance weighting
120 void interpolateInternalField
122 const GeometricField<Type, fvPatchField, volMesh>&,
123 GeometricField<Type, pointPatchField, pointMesh>&
126 //- Interpolate from volField to pointField
127 // using inverse distance weighting
131 const GeometricField<Type, fvPatchField, volMesh>&,
132 GeometricField<Type, pointPatchField, pointMesh>&
135 //- Interpolate volField using inverse distance weighting
136 // returning pointField with the same patchField types
138 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
140 const GeometricField<Type, fvPatchField, volMesh>&,
141 const wordList& patchFieldTypes
144 //- Interpolate tmp<volField> using inverse distance weighting
145 // returning pointField with the same patchField types
147 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
149 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
150 const wordList& patchFieldTypes
153 //- Interpolate volField using inverse distance weighting
154 // returning pointField
156 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
158 const GeometricField<Type, fvPatchField, volMesh>&
161 //- Interpolate tmp<volField> using inverse distance weighting
162 // returning pointField
164 tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
166 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 } // End namespace Foam
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 # include "volPointInterpolate.C"
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 // ************************************************************************* //