initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / interpolation / volPointInterpolation / volPointInterpolation.H
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
10
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.
15
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.
20
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
24
25 Class
26     Foam::volPointInterpolation
27
28 Description
29     Foam::volPointInterpolation
30
31 SourceFiles
32     volPointInterpolation.C
33     volPointInterpolate.C
34
35 \*---------------------------------------------------------------------------*/
36
37 #ifndef volPointInterpolation_H
38 #define volPointInterpolation_H
39
40 #include "MeshObject.H"
41 #include "pointPatchInterpolation.H"
42
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45 namespace Foam
46 {
47
48 class fvMesh;
49 class pointMesh;
50
51 /*---------------------------------------------------------------------------*\
52                        Class volPointInterpolation Declaration
53 \*---------------------------------------------------------------------------*/
54
55 class volPointInterpolation
56 :
57     public MeshObject<fvMesh, volPointInterpolation>
58 {
59     // Private data
60
61         //- Boundary interpolation engine.
62         pointPatchInterpolation boundaryInterpolator_;
63
64         //- Interpolation scheme weighting factor array.
65         scalarListList pointWeights_;
66
67
68     // Private member functions
69
70         //- Construct point weighting factors
71         void makeWeights();
72
73         //- Disallow default bitwise copy construct
74         volPointInterpolation(const volPointInterpolation&);
75
76         //- Disallow default bitwise assignment
77         void operator=(const volPointInterpolation&);
78
79
80 public:
81
82     // Declare name of the class and its debug switch
83     ClassName("volPointInterpolation");
84
85
86     // Constructors
87
88         //- Constructor given fvMesh and pointMesh.
89         explicit volPointInterpolation(const fvMesh&);
90
91
92     // Destructor
93
94         ~volPointInterpolation();
95
96
97     // Member functions
98
99         // Access
100
101             const fvMesh& mesh() const
102             {
103                 return boundaryInterpolator_.mesh();
104             }
105
106
107         // Edit
108
109             //- Update mesh topology using the morph engine
110             void updateMesh();
111
112             //- Correct weighting factors for moving mesh.
113             bool movePoints();
114
115
116     // Interpolation functions
117
118         //- Interpolate internal field from volField to pointField
119         //  using inverse distance weighting
120         template<class Type>
121         void interpolateInternalField
122         (
123             const GeometricField<Type, fvPatchField, volMesh>&,
124             GeometricField<Type, pointPatchField, pointMesh>&
125         ) const;
126
127         //- Interpolate from volField to pointField
128         //  using inverse distance weighting
129         template<class Type>
130         void interpolate
131         (
132             const GeometricField<Type, fvPatchField, volMesh>&,
133             GeometricField<Type, pointPatchField, pointMesh>&
134         ) const;
135
136         //- Interpolate volField using inverse distance weighting
137         //  returning pointField with the same patchField types
138         template<class Type>
139         tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
140         (
141             const GeometricField<Type, fvPatchField, volMesh>&,
142             const wordList& patchFieldTypes
143         ) const;
144
145         //- Interpolate tmp<volField> using inverse distance weighting
146         //  returning pointField with the same patchField types
147         template<class Type>
148         tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
149         (
150             const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
151             const wordList& patchFieldTypes
152         ) const;
153
154         //- Interpolate volField using inverse distance weighting 
155         //  returning pointField
156         template<class Type>
157         tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
158         (
159             const GeometricField<Type, fvPatchField, volMesh>&
160         ) const;
161
162         //- Interpolate tmp<volField> using inverse distance weighting 
163         //  returning pointField
164         template<class Type>
165         tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate
166         (
167             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
168         ) const;
169 };
170
171
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173
174 } // End namespace Foam
175
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177
178 #ifdef NoRepository
179 #   include "volPointInterpolate.C"
180 #endif
181
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183
184 #endif
185
186 // ************************************************************************* //