interpolation.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::interpolation
29 
30 Description
31  Abstract base class for volume field interpolation
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef Foam_interpolation_H
36 #define Foam_interpolation_H
37 
38 #include "faceList.H"
39 #include "volFieldsFwd.H"
40 #include "pointFields.H"
41 #include "autoPtr.H"
42 #include "typeInfo.H"
43 #include "barycentric.H"
44 #include "tetIndices.H"
45 #include "runTimeSelectionTables.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward Declarations
53 class polyMesh;
54 
55 /*---------------------------------------------------------------------------*\
56  Class interpolation Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class Type>
60 class interpolation
61 {
62 protected:
63 
64  // Protected Data
65 
67 
68  const polyMesh& pMesh_;
70  const faceList& pMeshFaces_;
73 
74 
75 public:
76 
77  //- Runtime type information
78  virtual const word& type() const = 0;
79 
80 
81  // Declare run-time constructor selection table
82 
84  (
85  autoPtr,
87  dictionary,
88  (
90  ),
91  (psi)
92  );
93 
94 
95  // Selectors
96 
97  //- Return a reference to the specified interpolation scheme
99  (
100  const word& interpolationType,
102  );
103 
104  //- Return a reference to the selected interpolation scheme
106  (
107  const dictionary& interpolationSchemes,
109  );
110 
111 
112  // Constructors
113 
114  //- Construct from components
115  explicit interpolation
116  (
118  );
119 
120 
121  //- Destructor
122  virtual ~interpolation() = default;
123 
124 
125  // Member Functions
126 
127  //- Return the field to be interpolated
129  {
130  return psi_;
131  }
132 
133  //- Interpolate field to the given point in the given cell
134  virtual Type interpolate
135  (
136  const vector& position,
137  const label celli,
138  const label facei = -1
139  ) const = 0;
140 
141  //- Interpolate field to the given coordinates in the tetrahedron
142  //- defined by the given indices.
143  // Calls interpolate function (vector, cell, face) except
144  // where overridden by derived interpolation types.
145  virtual Type interpolate
146  (
147  const barycentric& coordinates,
148  const tetIndices& tetIs,
149  const label facei = -1
150  ) const
151  {
152  return
154  (
156  tetIs.cell(),
157  facei
158  );
159  }
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #define makeInterpolationType(SS, Type) \
170  \
171 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
172  \
173 interpolation<Type>::adddictionaryConstructorToTable<SS<Type>> \
174  add##SS##Type##ConstructorToTable_;
175 
176 
177 #define makeInterpolation(SS) \
178  \
179 makeInterpolationType(SS, scalar) \
180 makeInterpolationType(SS, vector) \
181 makeInterpolationType(SS, sphericalTensor) \
182 makeInterpolationType(SS, symmTensor) \
183 makeInterpolationType(SS, tensor)
184 
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 #ifdef NoRepository
189  #include "interpolation.C"
190 #endif
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 #endif
195 
196 // ************************************************************************* //
const GeometricField< Type, fvPatchField, volMesh > & psi_
Definition: interpolation.H:59
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const polyMesh & pMesh_
Definition: interpolation.H:61
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
virtual const word & type() const =0
Runtime type information.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
Definition: tetIndicesI.H:129
Point barycentricToPoint(const barycentric &bary) const
Calculate the point from the given barycentric coordinates.
Definition: tetrahedronI.H:336
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 vectorField & pMeshFaceCentres_
Definition: interpolation.H:64
const vectorField & pMeshPoints_
Definition: interpolation.H:62
label cell() const noexcept
Return the cell index.
Definition: tetIndices.H:146
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual ~interpolation()=default
Destructor.
const faceList & pMeshFaces_
Definition: interpolation.H:63
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
declareRunTimeSelectionTable(autoPtr, interpolation, dictionary,(const GeometricField< Type, fvPatchField, volMesh > &psi),(psi))
PtrList< coordinateSystem > coordinates(solidRegions.size())
Abstract base class for volume field interpolation.
interpolation(const GeometricField< Type, fvPatchField, volMesh > &psi)
Construct from components.
Definition: interpolation.C:30
const vectorField & pMeshFaceAreas_
Definition: interpolation.H:65
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Macros to ease declaration of run-time selection tables.
Namespace for OpenFOAM.
const GeometricField< Type, fvPatchField, volMesh > & psi() const noexcept
Return the field to be interpolated.