31#include "demandDrivenData.H"
35const Foam::scalarSquareMatrix& Foam::RBFInterpolation::B()
const
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());
51 for (
int i = 0;
i <
A.n();
i++)
53 for (
int k = 0; k <
A.n(); k++)
60 Eigen::MatrixXd invEig = Aeig.fullPivHouseholderQr().solve(one);
62 for (
int i = 0;
i < invEig.rows();
i++)
64 for (
int k = 0; k < invEig.cols(); k++)
66 invMatrix[
i][k] = invEig(
i, k);
73void Foam::RBFInterpolation::calcB()
const
83 Eigen::MatrixXd Aeig = Eigen::MatrixXd::Zero(controlPoints_.size() + polySize,
84 controlPoints_.size() + polySize);
85 const label nControlPoints = controlPoints_.size();
87 for (label
i = 0;
i < nControlPoints;
i++)
89 scalarField weights = RBF_->weights(controlPoints_, controlPoints_[
i]);
91 for (label
col = 0;
col < nControlPoints;
col++)
101 label row = nControlPoints;
102 row < nControlPoints + 1;
106 for (label
col = 0;
col < nControlPoints;
col++)
108 Aeig(
col, row) = 1.0;
109 Aeig(row,
col) = 1.0;
116 label row = nControlPoints + 1;
117 row < nControlPoints + 2;
121 for (label
col = 0;
col < nControlPoints;
col++)
123 Aeig(
col, row) = controlPoints_[
col].x();
124 Aeig(row,
col) = controlPoints_[
col].x();
131 label row = nControlPoints + 2;
132 row < nControlPoints + 3;
136 for (label
col = 0;
col < nControlPoints;
col++)
138 Aeig(
col, row) = controlPoints_[
col].y();
139 Aeig(row,
col) = controlPoints_[
col].y();
146 label row = nControlPoints + 3;
147 row < nControlPoints + 4;
151 for (label
col = 0;
col < nControlPoints;
col++)
153 Aeig(
col, row) = controlPoints_[
col].z();
154 Aeig(row,
col) = controlPoints_[
col].z();
161 label row = nControlPoints;
162 row < nControlPoints + 4;
168 label
col = nControlPoints;
169 col < nControlPoints + 4;
173 Aeig(row,
col) = 0.0;
178 Info <<
"Inverting RBF motion matrix" << endl;
179 Eigen::MatrixXd InvAeig = Aeig.fullPivLu().inverse();
180 simpleMatrix<scalar> InvA(controlPoints_.size() + polySize);
182 for (
int i = 0;
i < InvAeig.rows();
i++)
184 for (
int k = 0; k < InvAeig.cols(); k++)
186 InvA[
i][k] = InvAeig(
i, k);
193 BPtr_ =
new scalarSquareMatrix(InvA);
197void Foam::RBFInterpolation::clearOut()
199 deleteDemandDrivenData(BPtr_);
207 const dictionary& dict,
208 const vectorField& controlPoints,
209 const vectorField& dataPoints
212 controlPoints_(controlPoints),
213 dataPoints_(dataPoints),
214 RBF_(
RBFFunction::New(word(dict.lookup(
"RBF")), dict)),
216 focalPoint_(dict.lookup(
"focalPoint")),
217 innerRadius_(readScalar(dict.lookup(
"innerRadius"))),
218 outerRadius_(readScalar(dict.lookup(
"outerRadius"))),
219 polynomials_(dict.lookup(
"polynomials"))
228 controlPoints_(rbf.controlPoints_),
229 dataPoints_(rbf.dataPoints_),
230 RBF_(rbf.RBF_->clone()),
232 focalPoint_(rbf.focalPoint_),
233 innerRadius_(rbf.innerRadius_),
234 outerRadius_(rbf.outerRadius_),
235 polynomials_(rbf.polynomials_)
RBFInterpolation(const dictionary &dict, const vectorField &controlPoints, const vectorField &dataPoints)
SquareMatrix< double > EigenInvert(SquareMatrix< double > &A)