libROM  v1.0
Data-driven physical simulation library
MatrixInterpolator.cpp
1 
11 // Description: Implementation of the MatrixInterpolator algorithm.
12 
13 #include "MatrixInterpolator.h"
14 
15 #include <limits.h>
16 #include <cmath>
17 #include "linalg/Matrix.h"
18 #include "linalg/scalapack_wrapper.h"
19 #include "mpi.h"
20 
21 /* Use C++11 built-in shared pointers if available; else fallback to Boost. */
22 #if __cplusplus >= 201103L
23 #include <memory>
24 #else
25 #include <boost/shared_ptr.hpp>
26 #endif
27 
28 /* Use automatically detected Fortran name-mangling scheme */
29 #define dposv CAROM_FC_GLOBAL(dposv, DPOSV)
30 
31 extern "C" {
32  // Solve a system of linear equations.
33  void dposv(char*, int*, int*, double*, int*, double*, int*, int*);
34 }
35 
36 using namespace std;
37 
38 namespace CAROM {
39 
40 MatrixInterpolator::MatrixInterpolator(const std::vector<Vector> &
41  parameter_points,
42  const std::vector<std::shared_ptr<Matrix>> & rotation_matrices,
43  const std::vector<std::shared_ptr<Matrix>> & reduced_matrices,
44  int ref_point,
45  std::string matrix_type,
46  std::string rbf,
47  std::string interp_method,
48  double closest_rbf_val,
49  bool compute_gradients) :
50  Interpolator(parameter_points,
51  rotation_matrices,
52  ref_point,
53  rbf,
54  interp_method,
55  closest_rbf_val,
56  compute_gradients)
57 {
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;
62 
63  // Rotate the reduced matrices
64  for (int i = 0; i < reduced_matrices.size(); i++)
65  {
66  // The ref_point does not need to be rotated
67  if (i == d_ref_point)
68  {
69  d_rotated_reduced_matrices.push_back(reduced_matrices[i]);
70  }
71  else
72  {
73  if (reduced_matrices[i]->numRows() == rotation_matrices[i]->numRows()
74  && reduced_matrices[i]->numColumns() == rotation_matrices[i]->numRows())
75  {
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);
80  }
81  else if (reduced_matrices[i]->numRows() == rotation_matrices[i]->numRows())
82  {
83  std::shared_ptr<Matrix> Q_tA = rotation_matrices[i]->transposeMult(
84  *reduced_matrices[i]);
85  d_rotated_reduced_matrices.push_back(Q_tA);
86  }
87  else if (reduced_matrices[i]->numColumns() == rotation_matrices[i]->numRows())
88  {
89  std::shared_ptr<Matrix> AQ = reduced_matrices[i]->mult(*rotation_matrices[i]);
90  d_rotated_reduced_matrices.push_back(AQ);
91  }
92  else
93  {
94  d_rotated_reduced_matrices.push_back(reduced_matrices[i]);
95  }
96  }
97 
98  }
99 }
100 
101 std::shared_ptr<Matrix> MatrixInterpolator::interpolate(const Vector & point,
102  bool orthogonalize)
103 {
104  std::shared_ptr<Matrix> interpolated_matrix;
105  if (d_matrix_type == "SPD")
106  {
108  {
109  CAROM_ERROR("Gradients are only implemented for B or G");
110  }
111  interpolated_matrix = interpolateSPDMatrix(point);
112  }
113  else if (d_matrix_type == "NS")
114  {
116  {
117  CAROM_ERROR("Gradients are only implemented for B or G");
118  }
119  interpolated_matrix = interpolateNonSingularMatrix(point);
120  }
121  else
122  {
123  interpolated_matrix = interpolateMatrix(point);
124  }
125 
126  // Orthogonalize the interpolated matrix if requested.
127  if (orthogonalize)
128  {
129  interpolated_matrix->orthogonalize();
130  }
131  return interpolated_matrix;
132 }
133 
134 void MatrixInterpolator::obtainLambda()
135 {
136  if (d_interp_method == "LS")
137  {
138  // Solving f = B*lambda
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++)
142  {
143  for (int j = 0; j < f_T->numColumns(); j++)
144  {
145  f_T->item(i, j) = d_gammas[j]->getData()[i];
146  }
147  }
148 
149  // Obtain B matrix by calculating RBF.
150  Matrix* B = new Matrix(d_gammas.size(), d_gammas.size(), false);
151  for (int i = 0; i < B->numRows(); i++)
152  {
153  B->item(i, i) = 1.0;
154  for (int j = i + 1; j < B->numColumns(); j++)
155  {
156  double res = obtainRBF(d_rbf, d_epsilon, d_parameter_points[i],
157  d_parameter_points[j]);
158  B->item(i, j) = res;
159  B->item(j, i) = res;
160  }
161  }
162 
163  char uplo = 'U';
164  int gamma_size = d_gammas.size();
165  int num_elements = d_gammas[0]->numRows() * d_gammas[0]->numColumns();
166  int info;
167 
168  dposv(&uplo, &gamma_size, &num_elements, B->getData(), &gamma_size,
169  f_T->getData(), &gamma_size, &info);
170  if (info != 0)
171  {
172  std::cout << "Linear solve failed. Please choose a different epsilon value." <<
173  std::endl;
174  }
175  CAROM_VERIFY(info == 0);
176 
177  delete B;
178 
179  d_lambda_T = std::unique_ptr<Matrix>(f_T);
180  }
181 }
182 
183 std::unique_ptr<Matrix> MatrixInterpolator::obtainLogInterpolatedMatrix(
184  std::vector<double>& rbf)
185 {
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());
190  if (d_interp_method == "LS")
191  {
192  for (int i = 0; i < d_lambda_T->numRows(); i++)
193  {
194  for (int j = 0; j < rbf.size(); j++)
195  {
196  log_interpolated_matrix->getData()[i] += d_lambda_T->item(i, j) * rbf[j];
197  }
198  }
199  }
200  else if (d_interp_method == "IDW")
201  {
202  double sum = rbfWeightedSum(rbf);
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++)
206  {
207  for (int j = 0; j < rbf.size(); j++)
208  {
209  log_interpolated_matrix->getData()[i] += d_gammas[j]->getData()[i] * rbf[j];
210  }
211  log_interpolated_matrix->getData()[i] /= sum;
212  }
213  }
214  else if (d_interp_method == "LP")
215  {
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++)
219  {
220  for (int j = 0; j < rbf.size(); j++)
221  {
222  log_interpolated_matrix->getData()[i] += d_gammas[j]->getData()[i] * rbf[j];
223  }
224  }
225  }
226  return std::unique_ptr<Matrix>(log_interpolated_matrix);
227 }
228 
229 std::shared_ptr<Matrix> MatrixInterpolator::interpolateSPDMatrix(
230  const Vector & point)
231 {
232  if (d_gammas.size() == 0)
233  {
234  // Diagonalize X to work towards obtaining X^-1/2
235  EigenPair ref_reduced_matrix_eigenpair = SymmetricRightEigenSolve(
236  *d_rotated_reduced_matrices[d_ref_point]);
237  Matrix* ref_reduced_matrix_ev = ref_reduced_matrix_eigenpair.ev;
238 
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++)
243  {
244  ref_reduced_matrix_sqrt_eigs->item(i, i) =
245  std::sqrt(ref_reduced_matrix_eigenpair.eigs[i]);
246  }
247 
248  Matrix ref_reduced_matrix_ev_inv(*ref_reduced_matrix_ev);
249  ref_reduced_matrix_ev->inverse(ref_reduced_matrix_ev_inv);
250 
251  // Obtain X^1/2
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);
256 
257  Matrix x_half_power_inv(*d_x_half_power);
258 
259  // Obtain X^-1/2
260  d_x_half_power->inverse(x_half_power_inv);
261 
262  delete ref_reduced_matrix_ev;
263  delete ref_reduced_matrix_sqrt_eigs;
264 
265  // Obtain gammas for all points in the database.
266  for (int i = 0; i < d_parameter_points.size(); i++)
267  {
268  // For the ref point, gamma is the zero matrix
269  if (i == d_ref_point)
270  {
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);
275  }
276  else
277  {
278  std::unique_ptr<Matrix> x_half_power_inv_mult_y = x_half_power_inv.mult(
279  *d_rotated_reduced_matrices[i]);
280 
281  // Obtain X^-1/2*Y*X^-1/2
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);
284 
285  // Diagonalize X^-1/2*Y*X^-1/2 to obtain the log of this matrix.
286  // Diagonalize YX^-1 to obtain log of this matrix.
287  // Following https://en.wikipedia.org/wiki/Logarithm_of_a_matrix
288  // 1. Diagonalize M to obtain M' = V^-1*M*V. M' are the eigenvalues
289  // of M and V are the eigenvectors of M.
290  // 2. log M = V(log M')V^-1
291  EigenPair log_eigenpair = SymmetricRightEigenSolve(
292  *x_half_power_inv_mult_y_mult_x_half_power_inv);
293  Matrix* log_ev = log_eigenpair.ev;
294 
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++)
298  {
299  if (log_eigenpair.eigs[i] < 0)
300  {
301  if (d_rank == 0) std::cout << "Some eigenvalues of this matrix are negative, "
302  <<
303  "which leads to NaN values when taking the log. Aborting." << std::endl;
304  CAROM_VERIFY(log_eigenpair.eigs[i] > 0);
305  }
306  log_eigs->item(i, i) = std::log(log_eigenpair.eigs[i]);
307  }
308 
309  // Invert matrix.
310  Matrix log_ev_inv(*log_ev);
311  log_ev->inverse(log_ev_inv);
312 
313  // Perform log mapping.
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);
317 
318  delete log_ev;
319  delete log_eigs;
320  }
321  }
322 
323  // Obtain lambda for the P interpolation matrix
324  obtainLambda();
325  }
326 
327  // Obtain distances from database points to new point
328  std::vector<double> rbf = obtainRBFToTrainingPoints(d_parameter_points,
330  d_rbf, d_epsilon, point);
331 
332  // Interpolate gammas to get gamma for new point
333  std::unique_ptr<Matrix> log_interpolated_matrix = obtainLogInterpolatedMatrix(
334  rbf);
335 
336  // Diagonalize the new gamma so we can exponentiate it
337  // Diagonalize X to obtain exp(X) of this matrix.
338  // Following https://en.wikipedia.org/wiki/Matrix_exponential
339  // 1. Diagonalize M to obtain M' = V^-1*M*V. M' are the eigenvalues
340  // of M and V are the eigenvectors of M.
341  // 2. exp M = V(exp M')V^-1
342  EigenPair exp_eigenpair = SymmetricRightEigenSolve(*log_interpolated_matrix);
343 
344  Matrix* exp_ev = exp_eigenpair.ev;
345 
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++)
349  {
350  exp_eigs->item(i, i) = std::exp(exp_eigenpair.eigs[i]);
351  }
352 
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);
356 
357  // Exponentiate gamma
358  std::unique_ptr<Matrix> exp_gamma = exp_ev_mult_exp_eig->mult(exp_ev_inv);
359 
360  delete exp_ev;
361  delete exp_eigs;
362 
363  // Obtain exp mapping by doing X^1/2*exp(gamma)*X^1/2
364  std::unique_ptr<Matrix> x_half_power_mult_exp_gamma = d_x_half_power->mult(
365  *exp_gamma);
366  std::shared_ptr<Matrix> interpolated_matrix = x_half_power_mult_exp_gamma->mult(
367  *d_x_half_power);
368 
369  return interpolated_matrix;
370 
371 }
372 
373 std::shared_ptr<Matrix> MatrixInterpolator::interpolateNonSingularMatrix(
374  const Vector & point)
375 {
376  if (d_gammas.size() == 0)
377  {
378  // Invert X
379  Matrix ref_matrix_inv(*d_rotated_reduced_matrices[d_ref_point]);
380  d_rotated_reduced_matrices[d_ref_point]->inverse(ref_matrix_inv);
381 
382  for (int i = 0; i < d_parameter_points.size(); i++)
383  {
384  // For ref_point, gamma is the zero matrix
385  if (i == d_ref_point)
386  {
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);
390  }
391  else
392  {
393  std::unique_ptr<Matrix> y_mult_ref_matrix_inv =
394  d_rotated_reduced_matrices[i]->mult(ref_matrix_inv);
395 
396  // Diagonalize YX^-1 to obtain log of this matrix.
397  // Following https://en.wikipedia.org/wiki/Logarithm_of_a_matrix
398  // 1. Diagonalize M to obtain M' = V^-1*M*V. M' are the eigenvalues
399  // of M and V are the eigenvectors of M.
400  // 2. log M = V(log M')V^-1
401  EigenPair log_eigenpair = SymmetricRightEigenSolve(*y_mult_ref_matrix_inv);
402  Matrix* log_ev = log_eigenpair.ev;
403 
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++)
407  {
408  if (log_eigenpair.eigs[i] < 0)
409  {
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." <<
412  std::endl;
413  CAROM_VERIFY(log_eigenpair.eigs[i] > 0);
414  }
415  log_eigs->item(i, i) = std::log(log_eigenpair.eigs[i]);
416  }
417 
418  Matrix log_ev_inv(log_ev->numRows(), log_ev->numColumns(), false);
419  log_ev->inverse(log_ev_inv);
420 
421  // Perform log mapping.
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);
425 
426  delete log_ev;
427  delete log_eigs;
428  }
429  }
430 
431  // Obtain lambda for the P interpolation matrix
432  obtainLambda();
433  }
434 
435  // Obtain distances from database points to new point
436  std::vector<double> rbf = obtainRBFToTrainingPoints(d_parameter_points,
437  d_interp_method, d_rbf, d_epsilon, point);
438 
439  // Interpolate gammas to get gamma for new point
440  std::unique_ptr<Matrix> log_interpolated_matrix = obtainLogInterpolatedMatrix(
441  rbf);
442 
443  // Diagonalize the new gamma so we can exponentiate it
444  // Diagonalize X to obtain exp(X) of this matrix.
445  // Following https://en.wikipedia.org/wiki/Matrix_exponential
446  // 1. Diagonalize M to obtain M' = V^-1*M*V. M' are the eigenvalues
447  // of M and V are the eigenvectors of M.
448  // 2. exp M = V(exp M')V^-1
449  EigenPair exp_eigenpair = SymmetricRightEigenSolve(*log_interpolated_matrix);
450  Matrix* exp_ev = exp_eigenpair.ev;
451  Matrix* exp_ev_inv = NULL;
452 
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++)
456  {
457  exp_eigs->item(i, i) = std::exp(exp_eigenpair.eigs[i]);
458  }
459 
460  // Invert matrix.
461  exp_ev->inverse(*exp_ev_inv);
462 
463  // Perform log mapping.
464  std::unique_ptr<Matrix> exp_ev_mult_exp_eig = exp_ev->mult(*exp_eigs);
465 
466  // Exponentiate gamma
467  std::unique_ptr<Matrix> exp_gamma = exp_ev_mult_exp_eig->mult(*exp_ev_inv);
468 
469  delete exp_ev;
470  delete exp_ev_inv;
471  delete exp_eigs;
472 
473  // Obtain exp mapping by doing exp(gamma)*X
474  std::shared_ptr<Matrix> interpolated_matrix = exp_gamma->mult(
475  *d_rotated_reduced_matrices[d_ref_point]);
476 
477  return interpolated_matrix;
478 }
479 
480 std::shared_ptr<Matrix> MatrixInterpolator::interpolateMatrix(
481  const Vector & point)
482 {
483  if (d_gammas.size() == 0)
484  {
485  for (int i = 0; i < d_parameter_points.size(); i++)
486  {
487  // For ref point, gamma is the zero matrix.
488  if (i == d_ref_point)
489  {
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);
495  }
496  else
497  {
498  // Gamma is Y - X
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);
502  }
503  }
504 
505  // Obtain lambda for the P interpolation matrix
506  obtainLambda();
507  }
508  // Obtain distances from database points to new point
509  std::vector<double> rbf = obtainRBFToTrainingPoints(d_parameter_points,
511  d_rbf, d_epsilon, point);
512 
513  // Interpolate gammas to get gamma for new point
514  std::shared_ptr<Matrix> interpolated_matrix(obtainLogInterpolatedMatrix(rbf));
515 
517  {
518  if(d_interp_method == "LS")
519  {
520  for (int i = 0; i < point.dim(); i++)
521  {
522  std::vector<double> rbf = obtainRBFGradientToTrainingPoints(d_parameter_points,
524  d_rbf, d_epsilon, point, i);
525 
526  std::shared_ptr<Matrix> gradient_matrix(obtainLogInterpolatedMatrix(rbf));
527  d_interpolation_gradient.push_back(gradient_matrix);
528  }
529  }
530  else
531  {
532  CAROM_ERROR("Interpolated gradients are only implemented for \"LS\"");
533  }
534  }
535 
536  // The exp mapping is X + the interpolated gamma
537  *interpolated_matrix += *d_rotated_reduced_matrices[d_ref_point];
538  return interpolated_matrix;
539 }
540 }
double d_epsilon
The RBF parameter that determines the width of influence. a small epsilon: larger influential width a...
Definition: Interpolator.h:95
int d_ref_point
The index within the vector of parameter points to the reference point.
Definition: Interpolator.h:76
int d_rank
The rank of the process this object belongs to.
Definition: Interpolator.h:65
std::unique_ptr< Matrix > d_lambda_T
The RHS of the linear solve in tangential space.
Definition: Interpolator.h:110
std::vector< Vector > d_parameter_points
The sampled parameter points.
Definition: Interpolator.h:100
bool d_compute_gradients
Flag that determines if a gradient with respect to the parameters should be computed.
Definition: Interpolator.h:116
std::string d_rbf
The RBF type (gaussian, multiquadric, inverse quadratic, inverse multiquadric)
Definition: Interpolator.h:82
std::string d_interp_method
The interpolation method (linear solve, inverse distance weighting, lagrangian polynomials)
Definition: Interpolator.h:88
double * getData() const
Get the matrix data as a pointer.
Definition: Matrix.h:844
int numColumns() const
Returns the number of columns in the Matrix. This method will return the same value from each process...
Definition: Matrix.h:226
const double & item(int row, int col) const
Const Matrix member access. Matrix data is stored in row-major format.
Definition: Matrix.h:746
int numRows() const
Returns the number of rows of the Matrix on this processor.
Definition: Matrix.h:200
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 > &parameter_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 > &parameter_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.
Definition: Matrix.cpp:1984