13 #include "MatrixInterpolator.h"
17 #include "linalg/Matrix.h"
18 #include "linalg/scalapack_wrapper.h"
22 #if __cplusplus >= 201103L
25 #include <boost/shared_ptr.hpp>
29 #define dposv CAROM_FC_GLOBAL(dposv, DPOSV)
33 void dposv(
char*,
int*,
int*,
double*,
int*,
double*,
int*,
int*);
40 MatrixInterpolator::MatrixInterpolator(
const std::vector<Vector> &
42 const std::vector<std::shared_ptr<Matrix>> & rotation_matrices,
43 const std::vector<std::shared_ptr<Matrix>> & reduced_matrices,
45 std::string matrix_type,
47 std::string interp_method,
48 double closest_rbf_val,
49 bool compute_gradients) :
58 CAROM_VERIFY(reduced_matrices.size() == rotation_matrices.size());
59 CAROM_VERIFY(matrix_type ==
"SPD" || matrix_type ==
"NS" || matrix_type ==
"R"
60 || matrix_type ==
"B");
61 d_matrix_type = matrix_type;
64 for (
int i = 0; i < reduced_matrices.size(); i++)
69 d_rotated_reduced_matrices.push_back(reduced_matrices[i]);
73 if (reduced_matrices[i]->numRows() == rotation_matrices[i]->numRows()
74 && reduced_matrices[i]->numColumns() == rotation_matrices[i]->numRows())
76 std::unique_ptr<Matrix> Q_tA = rotation_matrices[i]->transposeMult(
77 *reduced_matrices[i]);
78 std::shared_ptr<Matrix> Q_tAQ = Q_tA->mult(*rotation_matrices[i]);
79 d_rotated_reduced_matrices.push_back(Q_tAQ);
81 else if (reduced_matrices[i]->numRows() == rotation_matrices[i]->numRows())
83 std::shared_ptr<Matrix> Q_tA = rotation_matrices[i]->transposeMult(
84 *reduced_matrices[i]);
85 d_rotated_reduced_matrices.push_back(Q_tA);
87 else if (reduced_matrices[i]->numColumns() == rotation_matrices[i]->numRows())
89 std::shared_ptr<Matrix> AQ = reduced_matrices[i]->mult(*rotation_matrices[i]);
90 d_rotated_reduced_matrices.push_back(AQ);
94 d_rotated_reduced_matrices.push_back(reduced_matrices[i]);
104 std::shared_ptr<Matrix> interpolated_matrix;
105 if (d_matrix_type ==
"SPD")
109 CAROM_ERROR(
"Gradients are only implemented for B or G");
111 interpolated_matrix = interpolateSPDMatrix(point);
113 else if (d_matrix_type ==
"NS")
117 CAROM_ERROR(
"Gradients are only implemented for B or G");
119 interpolated_matrix = interpolateNonSingularMatrix(point);
123 interpolated_matrix = interpolateMatrix(point);
129 interpolated_matrix->orthogonalize();
131 return interpolated_matrix;
134 void MatrixInterpolator::obtainLambda()
139 Matrix* f_T =
new Matrix(d_gammas[0]->numRows() * d_gammas[0]->numColumns(),
140 d_gammas.size(),
false);
141 for (
int i = 0; i < f_T->
numRows(); i++)
145 f_T->
item(i, j) = d_gammas[j]->getData()[i];
150 Matrix* B =
new Matrix(d_gammas.size(), d_gammas.size(),
false);
151 for (
int i = 0; i < B->numRows(); i++)
154 for (
int j = i + 1; j < B->numColumns(); j++)
164 int gamma_size = d_gammas.size();
165 int num_elements = d_gammas[0]->numRows() * d_gammas[0]->numColumns();
168 dposv(&uplo, &gamma_size, &num_elements, B->getData(), &gamma_size,
169 f_T->
getData(), &gamma_size, &info);
172 std::cout <<
"Linear solve failed. Please choose a different epsilon value." <<
175 CAROM_VERIFY(info == 0);
183 std::unique_ptr<Matrix> MatrixInterpolator::obtainLogInterpolatedMatrix(
184 std::vector<double>& rbf)
186 Matrix* log_interpolated_matrix =
new Matrix(
187 d_rotated_reduced_matrices[
d_ref_point]->numRows(),
188 d_rotated_reduced_matrices[
d_ref_point]->numColumns(),
189 d_rotated_reduced_matrices[
d_ref_point]->distributed());
192 for (
int i = 0; i <
d_lambda_T->numRows(); i++)
194 for (
int j = 0; j < rbf.size(); j++)
196 log_interpolated_matrix->getData()[i] +=
d_lambda_T->item(i, j) * rbf[j];
203 int num_elements = d_rotated_reduced_matrices[
d_ref_point]->numRows() *
204 d_rotated_reduced_matrices[
d_ref_point]->numColumns();
205 for (
int i = 0; i < num_elements; i++)
207 for (
int j = 0; j < rbf.size(); j++)
209 log_interpolated_matrix->getData()[i] += d_gammas[j]->getData()[i] * rbf[j];
211 log_interpolated_matrix->getData()[i] /= sum;
216 int num_elements = d_rotated_reduced_matrices[
d_ref_point]->numRows() *
217 d_rotated_reduced_matrices[
d_ref_point]->numColumns();
218 for (
int i = 0; i < num_elements; i++)
220 for (
int j = 0; j < rbf.size(); j++)
222 log_interpolated_matrix->getData()[i] += d_gammas[j]->getData()[i] * rbf[j];
226 return std::unique_ptr<Matrix>(log_interpolated_matrix);
229 std::shared_ptr<Matrix> MatrixInterpolator::interpolateSPDMatrix(
230 const Vector & point)
232 if (d_gammas.size() == 0)
237 Matrix* ref_reduced_matrix_ev = ref_reduced_matrix_eigenpair.ev;
239 Matrix* ref_reduced_matrix_sqrt_eigs =
new Matrix(
240 ref_reduced_matrix_eigenpair.eigs.size(),
241 ref_reduced_matrix_eigenpair.eigs.size(),
false);
242 for (
int i = 0; i < ref_reduced_matrix_eigenpair.eigs.size(); i++)
244 ref_reduced_matrix_sqrt_eigs->item(i, i) =
245 std::sqrt(ref_reduced_matrix_eigenpair.eigs[i]);
248 Matrix ref_reduced_matrix_ev_inv(*ref_reduced_matrix_ev);
249 ref_reduced_matrix_ev->inverse(ref_reduced_matrix_ev_inv);
252 std::unique_ptr<Matrix> ref_reduced_matrix_ev_mult_sqrt_eig =
253 ref_reduced_matrix_ev->mult(*ref_reduced_matrix_sqrt_eigs);
254 d_x_half_power = ref_reduced_matrix_ev_mult_sqrt_eig->mult(
255 ref_reduced_matrix_ev_inv);
257 Matrix x_half_power_inv(*d_x_half_power);
260 d_x_half_power->inverse(x_half_power_inv);
262 delete ref_reduced_matrix_ev;
263 delete ref_reduced_matrix_sqrt_eigs;
271 std::shared_ptr<Matrix> gamma(
new Matrix(x_half_power_inv.numRows(),
272 x_half_power_inv.numColumns(),
273 x_half_power_inv.distributed()));
274 d_gammas.push_back(gamma);
278 std::unique_ptr<Matrix> x_half_power_inv_mult_y = x_half_power_inv.mult(
279 *d_rotated_reduced_matrices[i]);
282 std::unique_ptr<Matrix> x_half_power_inv_mult_y_mult_x_half_power_inv =
283 x_half_power_inv_mult_y->mult(x_half_power_inv);
292 *x_half_power_inv_mult_y_mult_x_half_power_inv);
293 Matrix* log_ev = log_eigenpair.ev;
295 Matrix* log_eigs =
new Matrix(log_eigenpair.eigs.size(),
296 log_eigenpair.eigs.size(),
false);
297 for (
int i = 0; i < log_eigenpair.eigs.size(); i++)
299 if (log_eigenpair.eigs[i] < 0)
301 if (
d_rank == 0) std::cout <<
"Some eigenvalues of this matrix are negative, "
303 "which leads to NaN values when taking the log. Aborting." << std::endl;
304 CAROM_VERIFY(log_eigenpair.eigs[i] > 0);
306 log_eigs->item(i, i) = std::log(log_eigenpair.eigs[i]);
310 Matrix log_ev_inv(*log_ev);
311 log_ev->inverse(log_ev_inv);
314 std::unique_ptr<Matrix> log_ev_mult_log_eig = log_ev->mult(*log_eigs);
315 std::shared_ptr<Matrix> gamma = log_ev_mult_log_eig->mult(log_ev_inv);
316 d_gammas.push_back(gamma);
333 std::unique_ptr<Matrix> log_interpolated_matrix = obtainLogInterpolatedMatrix(
344 Matrix* exp_ev = exp_eigenpair.ev;
346 Matrix* exp_eigs =
new Matrix(exp_eigenpair.eigs.size(),
347 exp_eigenpair.eigs.size(),
false);
348 for (
int i = 0; i < exp_eigenpair.eigs.size(); i++)
350 exp_eigs->item(i, i) = std::exp(exp_eigenpair.eigs[i]);
353 Matrix exp_ev_inv(*exp_ev);
354 exp_ev->inverse(exp_ev_inv);
355 std::unique_ptr<Matrix> exp_ev_mult_exp_eig = exp_ev->mult(*exp_eigs);
358 std::unique_ptr<Matrix> exp_gamma = exp_ev_mult_exp_eig->mult(exp_ev_inv);
364 std::unique_ptr<Matrix> x_half_power_mult_exp_gamma = d_x_half_power->mult(
366 std::shared_ptr<Matrix> interpolated_matrix = x_half_power_mult_exp_gamma->mult(
369 return interpolated_matrix;
373 std::shared_ptr<Matrix> MatrixInterpolator::interpolateNonSingularMatrix(
374 const Vector & point)
376 if (d_gammas.size() == 0)
379 Matrix ref_matrix_inv(*d_rotated_reduced_matrices[
d_ref_point]);
380 d_rotated_reduced_matrices[
d_ref_point]->inverse(ref_matrix_inv);
387 std::shared_ptr<Matrix> gamma(
new Matrix(ref_matrix_inv.numRows(),
388 ref_matrix_inv.numColumns(), ref_matrix_inv.distributed()));
389 d_gammas.push_back(gamma);
393 std::unique_ptr<Matrix> y_mult_ref_matrix_inv =
394 d_rotated_reduced_matrices[i]->mult(ref_matrix_inv);
402 Matrix* log_ev = log_eigenpair.ev;
404 Matrix* log_eigs =
new Matrix(log_eigenpair.eigs.size(),
405 log_eigenpair.eigs.size(),
false);
406 for (
int i = 0; i < log_eigenpair.eigs.size(); i++)
408 if (log_eigenpair.eigs[i] < 0)
410 if (
d_rank == 0) std::cout <<
"Some eigenvalues of this matrix are " <<
411 "negative, which leads to NaN values when taking the log. Aborting." <<
413 CAROM_VERIFY(log_eigenpair.eigs[i] > 0);
415 log_eigs->item(i, i) = std::log(log_eigenpair.eigs[i]);
418 Matrix log_ev_inv(log_ev->numRows(), log_ev->numColumns(),
false);
419 log_ev->inverse(log_ev_inv);
422 std::unique_ptr<Matrix> log_ev_mult_log_eig = log_ev->mult(*log_eigs);
423 std::shared_ptr<Matrix> gamma = log_ev_mult_log_eig->mult(log_ev_inv);
424 d_gammas.push_back(gamma);
440 std::unique_ptr<Matrix> log_interpolated_matrix = obtainLogInterpolatedMatrix(
450 Matrix* exp_ev = exp_eigenpair.ev;
451 Matrix* exp_ev_inv = NULL;
453 Matrix* exp_eigs =
new Matrix(exp_eigenpair.eigs.size(),
454 exp_eigenpair.eigs.size(),
false);
455 for (
int i = 0; i < exp_eigenpair.eigs.size(); i++)
457 exp_eigs->item(i, i) = std::exp(exp_eigenpair.eigs[i]);
461 exp_ev->inverse(*exp_ev_inv);
464 std::unique_ptr<Matrix> exp_ev_mult_exp_eig = exp_ev->mult(*exp_eigs);
467 std::unique_ptr<Matrix> exp_gamma = exp_ev_mult_exp_eig->mult(*exp_ev_inv);
474 std::shared_ptr<Matrix> interpolated_matrix = exp_gamma->mult(
477 return interpolated_matrix;
480 std::shared_ptr<Matrix> MatrixInterpolator::interpolateMatrix(
481 const Vector & point)
483 if (d_gammas.size() == 0)
490 std::shared_ptr<Matrix> gamma(
new Matrix(
491 d_rotated_reduced_matrices[
d_ref_point]->numRows(),
492 d_rotated_reduced_matrices[
d_ref_point]->numColumns(),
493 d_rotated_reduced_matrices[
d_ref_point]->distributed()));
494 d_gammas.push_back(gamma);
499 std::shared_ptr<Matrix> gamma(
new Matrix(*d_rotated_reduced_matrices[i]));
500 *gamma -= *d_rotated_reduced_matrices[
d_ref_point];
501 d_gammas.push_back(gamma);
514 std::shared_ptr<Matrix> interpolated_matrix(obtainLogInterpolatedMatrix(rbf));
520 for (
int i = 0; i < point.dim(); i++)
526 std::shared_ptr<Matrix> gradient_matrix(obtainLogInterpolatedMatrix(rbf));
527 d_interpolation_gradient.push_back(gradient_matrix);
532 CAROM_ERROR(
"Interpolated gradients are only implemented for \"LS\"");
537 *interpolated_matrix += *d_rotated_reduced_matrices[
d_ref_point];
538 return interpolated_matrix;
double d_epsilon
The RBF parameter that determines the width of influence. a small epsilon: larger influential width a...
int d_ref_point
The index within the vector of parameter points to the reference point.
int d_rank
The rank of the process this object belongs to.
std::unique_ptr< Matrix > d_lambda_T
The RHS of the linear solve in tangential space.
std::vector< Vector > d_parameter_points
The sampled parameter points.
bool d_compute_gradients
Flag that determines if a gradient with respect to the parameters should be computed.
std::string d_rbf
The RBF type (gaussian, multiquadric, inverse quadratic, inverse multiquadric)
std::string d_interp_method
The interpolation method (linear solve, inverse distance weighting, lagrangian polynomials)
double * getData() const
Get the matrix data as a pointer.
int numColumns() const
Returns the number of columns in the Matrix. This method will return the same value from each process...
const double & item(int row, int col) const
Const Matrix member access. Matrix data is stored in row-major format.
int numRows() const
Returns the number of rows of the Matrix on this processor.
std::shared_ptr< Matrix > interpolate(const Vector &point, bool orthogonalize=false)
Obtain the interpolated reduced matrix of the unsampled parameter point.
double rbfWeightedSum(std::vector< double > &rbf)
Compute the sum of the RBF weights.
std::vector< double > obtainRBFGradientToTrainingPoints(const std::vector< Vector > ¶meter_points, const std::string &interp_method, const std::string &rbf, double epsilon, const Vector &point, const int index)
Compute the RBF gradient from the parameter points with the unsampled parameter point.
std::vector< double > obtainRBFToTrainingPoints(const std::vector< Vector > ¶meter_points, const std::string &interp_method, const std::string &rbf, double epsilon, const Vector &point)
Compute the RBF from the parameter points with the unsampled parameter point.
double obtainRBF(std::string rbf, double epsilon, const Vector &point1, const Vector &point2)
Compute the RBF between two points.
EigenPair SymmetricRightEigenSolve(const Matrix &A)
Computes the eigenvectors/eigenvalues of an NxN real symmetric matrix.