interpolation.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your 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, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::interpolation
26 
27 Description
28  Abstract base class for interpolation
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef interpolation_H
33 #define interpolation_H
34 
35 #include "faceList.H"
36 #include "volFieldsFwd.H"
37 #include "pointFields.H"
38 #include "typeInfo.H"
39 #include "autoPtr.H"
40 #include "runTimeSelectionTables.H"
41 #include "tetIndices.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 class polyMesh;
49 
50 /*---------------------------------------------------------------------------*\
51  Class interpolation Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 template<class Type>
55 class interpolation
56 {
57 protected:
58 
59  // Protected data
60 
61  //- The vol field to interpolate
63 
64  //- Reference to the mesh
65  const polyMesh& mesh_;
66 
67 
68 public:
69 
70  //- Runtime type information
71  virtual const word& type() const = 0;
72 
73 
74  // Declare run-time constructor selection table
75 
77  (
78  autoPtr,
80  dictionary,
81  (
83  ),
84  (psi)
85  );
86 
87 
88  // Selectors
89 
90  //- Return a reference to the specified interpolation scheme
92  (
93  const word& interpolationType,
95  );
96 
97  //- Return a reference to the selected interpolation scheme
99  (
100  const dictionary& interpolationSchemes,
102  );
103 
104 
105  // Constructors
106 
107  //- Construct from components
109  (
111  );
112 
113 
114  //- Destructor
115  virtual ~interpolation()
116  {}
117 
118 
119  // Member Functions
120 
121  //- Return the field to be interpolated
123  {
124  return psi_;
125  }
126 
127  //- Interpolate field to the given point in the given cell
128  virtual Type interpolate
129  (
130  const vector& position,
131  const label celli,
132  const label facei = -1
133  ) const = 0;
134 
135  //- As above, but for a field
137  (
138  const vectorField& position,
139  const labelField& celli,
140  const labelField& facei = NullObjectRef<labelField>()
141  ) const = 0;
142 
143  //- Interpolate field to the given coordinates in the tetrahedron
144  // defined by the given indices. Calls interpolate function
145  // above here except where overridden by derived
146  // interpolation types.
147  virtual Type interpolate
148  (
149  const barycentric& coordinates,
150  const tetIndices& tetIs,
151  const label facei = -1
152  ) const;
153 
154  //- As above, but for a field
156  (
157  const Field<barycentric>& coordinates,
158  const labelField& celli,
159  const labelField& tetFacei,
160  const labelField& tetPti,
161  const labelField& facei = NullObjectRef<labelField>()
162  ) const = 0;
163 };
164 
165 
166 /*---------------------------------------------------------------------------*\
167  Class fieldInterpolation Declaration
168 \*---------------------------------------------------------------------------*/
169 
170 template<class Type, class InterpolationType>
171 class fieldInterpolation
172 :
173  public interpolation<Type>
174 {
175 public:
176 
177  // Constructors
178 
179  //- Inherit base class constructors
181 
182 
183  // Member Functions
184 
185  //- Interpolate field to the given points in the given cells
187  (
188  const vectorField& position,
189  const labelField& celli,
190  const labelField& facei = NullObjectRef<labelField>()
191  ) const;
192 
193  //- Interpolate field to the given coordinates in the given tetrahedra
195  (
196  const Field<barycentric>& coordinates,
197  const labelField& celli,
198  const labelField& tetFacei,
199  const labelField& tetPti,
200  const labelField& facei = NullObjectRef<labelField>()
201  ) const;
202 };
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace Foam
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 #define makeInterpolationType(SS, Type) \
212  \
213 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
214  \
215 interpolation<Type>::adddictionaryConstructorToTable<SS<Type>> \
216  add##SS##Type##ConstructorToTable_;
217 
219 #define makeInterpolation(SS) \
220  \
221 makeInterpolationType(SS, scalar) \
222 makeInterpolationType(SS, vector) \
223 makeInterpolationType(SS, sphericalTensor) \
224 makeInterpolationType(SS, symmTensor) \
225 makeInterpolationType(SS, tensor)
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 #ifdef NoRepository
231  #include "interpolation.C"
232 #endif
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 #endif
237 
238 // ************************************************************************* //
const GeometricField< Type, fvPatchField, volMesh > & psi_
The vol field to interpolate.
Definition: interpolation.H:61
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolation.C:48
virtual ~interpolation()
Destructor.
virtual const word & type() const =0
Runtime type information.
Generic GeometricField class.
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Return the field to be interpolated.
A class for handling words, derived from string.
Definition: word.H:59
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
const polyMesh & mesh_
Reference to the mesh.
Definition: interpolation.H:64
declareRunTimeSelectionTable(autoPtr, interpolation, dictionary,(const GeometricField< Type, fvPatchField, volMesh > &psi),(psi))
Abstract base class for interpolation.
interpolation(const GeometricField< Type, fvPatchField, volMesh > &psi)
Construct from components.
Definition: interpolation.C:35
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.