Loading...
Searching...
No Matches
RBFInterpolation.C
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
24Author
25 Frank Bos, TU Delft. All rights reserved.
26 Dubravko Matijasevic, FSB Zagreb.
27
28\*---------------------------------------------------------------------------*/
29
30#include "RBFInterpolation.H"
31#include "demandDrivenData.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35const Foam::scalarSquareMatrix& Foam::RBFInterpolation::B() const
36{
37 if (!BPtr_)
38 {
39 calcB();
40 }
41
42 return *BPtr_;
43}
44
45Foam::SquareMatrix<double> Foam::EigenInvert(Foam::SquareMatrix<double>& A)
46{
47 Foam::SquareMatrix<double> invMatrix = A;
48 Eigen::MatrixXd Aeig(A.n(), A.n());
49 Eigen::MatrixXd one = Eigen::MatrixXd::Identity(A.n(), A.n());
50
51 for (int i = 0; i < A.n(); i++)
52 {
53 for (int k = 0; k < A.n(); k++)
54 {
55 Aeig(i, k) = A[i][k];
56 }
57 }
58
59 //Eigen::MatrixXd invEig = Aeig.inverse();
60 Eigen::MatrixXd invEig = Aeig.fullPivHouseholderQr().solve(one);
61
62 for (int i = 0; i < invEig.rows(); i++)
63 {
64 for (int k = 0; k < invEig.cols(); k++)
65 {
66 invMatrix[i][k] = invEig(i, k);
67 }
68 }
69
70 return invMatrix;
71}
72
73void Foam::RBFInterpolation::calcB() const
74{
75 // Determine inverse of boundary connectivity matrix
76 label polySize(4);
77
78 if (!polynomials_)
79 {
80 polySize = 0;
81 }
82
83 Eigen::MatrixXd Aeig = Eigen::MatrixXd::Zero(controlPoints_.size() + polySize,
84 controlPoints_.size() + polySize);
85 const label nControlPoints = controlPoints_.size();
86
87 for (label i = 0; i < nControlPoints; i++)
88 {
89 scalarField weights = RBF_->weights(controlPoints_, controlPoints_[i]);
90
91 for (label col = 0; col < nControlPoints; col++)
92 {
93 Aeig(i, col) = weights[col];
94 }
95 }
96
97 if (polynomials_)
98 {
99 for
100 (
101 label row = nControlPoints;
102 row < nControlPoints + 1;
103 row++
104 )
105 {
106 for (label col = 0; col < nControlPoints; col++)
107 {
108 Aeig(col, row) = 1.0;
109 Aeig(row, col) = 1.0;
110 }
111 }
112
113 // Fill in X components of polynomial part of matrix
114 for
115 (
116 label row = nControlPoints + 1;
117 row < nControlPoints + 2;
118 row++
119 )
120 {
121 for (label col = 0; col < nControlPoints; col++)
122 {
123 Aeig(col, row) = controlPoints_[col].x();
124 Aeig(row, col) = controlPoints_[col].x();
125 }
126 }
127
128 // Fill in Y components of polynomial part of matrix
129 for
130 (
131 label row = nControlPoints + 2;
132 row < nControlPoints + 3;
133 row++
134 )
135 {
136 for (label col = 0; col < nControlPoints; col++)
137 {
138 Aeig(col, row) = controlPoints_[col].y();
139 Aeig(row, col) = controlPoints_[col].y();
140 }
141 }
142
143 // Fill in Z components of polynomial part of matrix
144 for
145 (
146 label row = nControlPoints + 3;
147 row < nControlPoints + 4;
148 row++
149 )
150 {
151 for (label col = 0; col < nControlPoints; col++)
152 {
153 Aeig(col, row) = controlPoints_[col].z();
154 Aeig(row, col) = controlPoints_[col].z();
155 }
156 }
157
158 // Fill 4x4 zero part of matrix
159 for
160 (
161 label row = nControlPoints;
162 row < nControlPoints + 4;
163 row++
164 )
165 {
166 for
167 (
168 label col = nControlPoints;
169 col < nControlPoints + 4;
170 col++
171 )
172 {
173 Aeig(row, col) = 0.0;
174 }
175 }
176 }
177
178 Info << "Inverting RBF motion matrix" << endl;
179 Eigen::MatrixXd InvAeig = Aeig.fullPivLu().inverse();
180 simpleMatrix<scalar> InvA(controlPoints_.size() + polySize);
181
182 for (int i = 0; i < InvAeig.rows(); i++)
183 {
184 for (int k = 0; k < InvAeig.cols(); k++)
185 {
186 InvA[i][k] = InvAeig(i, k);
187 }
188 }
189
190 // HJ and FB (05 Jan 2009)
191 // Collect ALL control points from ALL CPUs
192 // Create an identical inverse for all CPUs
193 BPtr_ = new scalarSquareMatrix(InvA);
194}
195
196
197void Foam::RBFInterpolation::clearOut()
198{
199 deleteDemandDrivenData(BPtr_);
200}
201
202
203// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
204
206(
207 const dictionary& dict,
208 const vectorField& controlPoints,
209 const vectorField& dataPoints
210)
211 :
212 controlPoints_(controlPoints),
213 dataPoints_(dataPoints),
214 RBF_(RBFFunction::New(word(dict.lookup("RBF")), dict)),
215 BPtr_(NULL),
216 focalPoint_(dict.lookup("focalPoint")),
217 innerRadius_(readScalar(dict.lookup("innerRadius"))),
218 outerRadius_(readScalar(dict.lookup("outerRadius"))),
219 polynomials_(dict.lookup("polynomials"))
220{}
221
222
224(
225 const RBFInterpolation& rbf
226)
227 :
228 controlPoints_(rbf.controlPoints_),
229 dataPoints_(rbf.dataPoints_),
230 RBF_(rbf.RBF_->clone()),
231 BPtr_(NULL),
232 focalPoint_(rbf.focalPoint_),
233 innerRadius_(rbf.innerRadius_),
234 outerRadius_(rbf.outerRadius_),
235 polynomials_(rbf.polynomials_)
236{}
237
238
239// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
240
245
246
247// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
248
250{
251 clearOut();
252}
253
254
255// ************************************************************************* //
RBFInterpolation(const dictionary &dict, const vectorField &controlPoints, const vectorField &dataPoints)
volScalarField & A
SquareMatrix< double > EigenInvert(SquareMatrix< double > &A)
singVal col(i)
label i
Definition pEqn.H:46