Loading...
Searching...
No Matches
RBFInterpolation.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 4.0
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7-------------------------------------------------------------------------------
8License
9 This file is part of foam-extend.
10
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.
15
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.
20
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/>.
23
24Class
25 RBFInterpolation
26
27Description
28 Radial basis function interpolation class
29
30Description
31 Interpolation class which uses Radial Basis Functions to interpolate
32 field from given data points to arbitrary set of points.
33
34 The coefficient vectors, alpha and beta are determined by solving
35 the system:
36
37 | db | = | Mbb Pb | | alpha |
38 | 0 | | Pb 0 | | beta |
39
40 where db are the given field values at data carrier points.
41 Mbb the boundary RBF correlation matrix (NbxNb), containing RBF evaluations
42 at the boundary nodes, and Pb some linear polynomial matrix (Nbx4).
43
44 In cases where far field data is not of interest, a cutoff function
45 is used to eliminate unnecessary data points in the far field
46
47Author
48 Frank Bos, TU Delft. All rights reserved.
49 Dubravko Matijasevic, FSB Zagreb.
50 Reorganisation by Hrvoje Jasak, Wikki Ltd.
51
52SourceFiles
53 RBFInterpolation.C
54 RBFInterpolationTemplates.C
55
56\*---------------------------------------------------------------------------*/
57
58#ifndef RBFInterpolation_H
59#define RBFInterpolation_H
60
61#include "dictionary.H"
62#include "RBFFunction.H"
63#include "point.H"
64#include "Switch.H"
65#include "simpleMatrix.H"
66#pragma GCC diagnostic push
67#pragma GCC diagnostic ignored "-Wold-style-cast"
68#pragma GCC diagnostic ignored "-Wnon-virtual-dtor"
69#include <Eigen/Eigen>
70#pragma GCC diagnostic pop
71// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72
73namespace Foam
74{
75
76/*---------------------------------------------------------------------------*\
77 Class RBFInterpolation Declaration
78\*---------------------------------------------------------------------------*/
79
80SquareMatrix<double> EigenInvert(SquareMatrix<double>& A);
81
83{
84 // Private data
85
86 //- Reference to control points
87 const vectorField& controlPoints_;
88
89 //- Reference to all points
90 const vectorField& dataPoints_;
91
92 //- RBF function
93 autoPtr<RBFFunction> RBF_;
94
95 //- Interpolation matrix
96 mutable scalarSquareMatrix* BPtr_;
97
98 //- Focal point for cut-off radii
99 point focalPoint_;
100
101 //- Inner cut-off radius
102 scalar innerRadius_;
103
104 //- Outer cut-off radius
105 scalar outerRadius_;
106
107 //- Add polynomials to RBF matrix
108 Switch polynomials_;
109
110
111 // Private Member Functions
112
113 //- Disallow default bitwise assignment
114 void operator=(const RBFInterpolation&);
115
116
117 //- Return interpolation matrix
118 const scalarSquareMatrix& B() const;
119
120 //- Calculate interpolation matrix
121 void calcB() const;
122
123 //- Clear out
124 void clearOut();
125
126
127 public:
128
129 // Constructors
130
131 //- Construct from components
133 (
134 const dictionary& dict,
135 const vectorField& controlPoints,
136 const vectorField& dataPoints
137 );
138
139 //- Construct as copy
141
142
143 // Destructor
144
146
147
148 // Member Functions
149
150 //- Return reference to control points
151 const vectorField& controlPoints() const
152 {
153 return controlPoints_;
154 }
155
156 //- Reference to all points
157 const vectorField& dataPoints() const
158 {
159 return dataPoints_;
160 }
161
162 //- Interpolate
163 template<class Type>
164 tmp<Field<Type>> interpolate(const Field<Type>& ctrlField) const;
165
166 //- Move points
167 void movePoints();
168};
169
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173} // End namespace Foam
174
175// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176
177#ifdef NoRepository
179#endif
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183#endif
184
185// ************************************************************************* //
tmp< Field< Type > > interpolate(const Field< Type > &ctrlField) const
RBFInterpolation(const dictionary &dict, const vectorField &controlPoints, const vectorField &dataPoints)
const vectorField & dataPoints() const
const vectorField & controlPoints() const
volScalarField & A
SquareMatrix< double > EigenInvert(SquareMatrix< double > &A)