libROM  v1.0
Data-driven physical simulation library
S_OPT.cpp
1 
11 // Description: Interface to the S_OPT algorithm to determine the rows of the
12 // basis to be sampled for the interpolation of the basis.
13 
14 #include "linalg/Matrix.h"
15 #include "Utilities.h"
16 #include "mpi.h"
17 #include <cmath>
18 #include <vector>
19 #include <map>
20 #include <set>
21 #include <cfloat>
22 
23 #include "S_OPT.h"
24 
25 using namespace std;
26 
27 namespace CAROM {
28 
29 void
30 S_OPT(const Matrix & f_basis,
31  int num_f_basis_vectors_used,
32  std::vector<int>& f_sampled_row,
33  std::vector<int>& f_sampled_rows_per_proc,
34  Matrix& f_basis_sampled_inv,
35  const int myid,
36  const int num_procs,
37  const int num_samples_req,
38  std::vector<int>* init_samples,
39  bool qr_factorize)
40 {
41  CAROM_VERIFY(num_procs == f_sampled_rows_per_proc.size());
42  // This algorithm determines the rows of f that should be sampled, the
43  // processor that owns each sampled row, and fills f_basis_sampled_inv with
44  // the inverse of the sampled rows of the basis.
45 
46  // Create an MPI_Datatype for the RowInfo struct.
47  MPI_Datatype MaxRowType, oldtypes[2];
48  int blockcounts[2];
49  MPI_Aint offsets[2], extent, lower_bound;
50  MPI_Status stat;
51  offsets[0] = 0;
52  oldtypes[0] = MPI_DOUBLE;
53  blockcounts[0] = 1;
54 
55  /*
56  MPI_Type_extent is deprecated in MPI-2 and removed in
57  MPI-3. MPI_Type_extent has been replaced by MPI_Type_get_extent.
58  */
59 #if (MPI_VERSION == 1)
60  MPI_Type_extent(MPI_DOUBLE, &extent);
61 #else
62  MPI_Type_get_extent(MPI_DOUBLE, &lower_bound, &extent);
63 #endif
64 
65  offsets[1] = extent;
66  oldtypes[1] = MPI_INT;
67  blockcounts[1] = 2;
68 
69  /*
70  MPI_Type_struct is deprecated in MPI-2 and removed in
71  MPI-3. MPI_Type_struct has been replaced by MPI_Type_create_struct.
72  */
73 #if (MPI_VERSION == 1)
74  MPI_Type_struct(2, blockcounts, offsets, oldtypes, &MaxRowType);
75 #else
76  MPI_Type_create_struct(2, blockcounts, offsets, oldtypes, &MaxRowType);
77 #endif
78  MPI_Type_commit(&MaxRowType);
79 
80  // Create an MPI_Op for the RowInfoMax function.
81  MPI_Op RowInfoOp;
82  MPI_Op_create((MPI_User_function*)RowInfoMax, true, &RowInfoOp);
83 
84  // Get the number of basis vectors and the size of each basis vector.
85  CAROM_VERIFY(0 < num_f_basis_vectors_used
86  && num_f_basis_vectors_used <= f_basis.numColumns());
87  const int num_basis_vectors =
88  std::min(num_f_basis_vectors_used, f_basis.numColumns());
89  const int num_samples = num_samples_req > 0 ? num_samples_req :
90  num_basis_vectors;
91  CAROM_VERIFY(num_basis_vectors <= num_samples
92  && num_samples <= f_basis.numDistributedRows());
93  CAROM_VERIFY(num_samples == f_sampled_row.size());
94  CAROM_VERIFY(num_samples == f_basis_sampled_inv.numRows() &&
95  num_basis_vectors == f_basis_sampled_inv.numColumns());
96  CAROM_VERIFY(!f_basis_sampled_inv.distributed());
97 
98  int num_rows = f_basis.numRows();
99 
100  // If num_basis_vectors is less than the number of columns of the basis,
101  // we need to truncate the basis.
102  std::shared_ptr<Matrix> f_basis_truncated_N;
103  const Matrix *f_basis_truncated = &f_basis;
104  if (num_basis_vectors < f_basis.numColumns())
105  {
106  f_basis_truncated_N = f_basis.getFirstNColumns(num_basis_vectors);
107  f_basis_truncated = f_basis_truncated_N.get();
108  }
109 
110  std::shared_ptr<Matrix> Vo_qr;
111  const Matrix *Vo = f_basis_truncated;
112  // Use the QR factorization of the input matrix, if requested
113  if (qr_factorize)
114  {
115  Vo_qr = f_basis_truncated->qr_factorize();
116  Vo = Vo_qr.get();
117  }
118 
119  int num_samples_obtained = 0;
120 
121  // Scratch space used throughout the algorithm.
122  std::vector<double> c(num_basis_vectors);
123 
124  vector<set<int> > proc_sampled_f_row(num_procs);
125  vector<map<int, int> > proc_f_row_to_tmp_fs_row(num_procs);
126  int num_f_basis_cols = f_basis_sampled_inv.numColumns();
127  Matrix V1(num_samples,
128  num_basis_vectors,
129  false);
130 
131  RowInfo f_bv_max_local, f_bv_max_global;
132 
133  // Gather information about initial samples given as input.
134  const int num_init_samples = init_samples ? init_samples->size() : 0;
135  int total_num_init_samples = 0;
136  MPI_Allreduce(&num_init_samples, &total_num_init_samples, 1,
137  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
138  CAROM_VERIFY(num_samples >= total_num_init_samples);
139 
140  int init_sample_offset = 0;
141  for (int i = 0; i < total_num_init_samples; i++)
142  {
143  f_bv_max_local.row_val = -DBL_MAX;
144  f_bv_max_local.proc = myid;
145  if (init_sample_offset < num_init_samples)
146  {
147  f_bv_max_local.row_val = 1.0;
148  f_bv_max_local.row = (*init_samples)[init_sample_offset];
149  }
150  MPI_Allreduce(&f_bv_max_local, &f_bv_max_global, 1,
151  MaxRowType, RowInfoOp, MPI_COMM_WORLD);
152  // Now get the first sampled row of the basis
153  if (f_bv_max_global.proc == myid) {
154  for (int j = 0; j < num_basis_vectors; ++j) {
155  c[j] = Vo->item(f_bv_max_global.row, j);
156  }
157  init_sample_offset++;
158  }
159  MPI_Bcast(c.data(), num_basis_vectors, MPI_DOUBLE,
160  f_bv_max_global.proc, MPI_COMM_WORLD);
161  // Now add the first sampled row of the basis to tmp_fs.
162  for (int j = 0; j < num_basis_vectors; ++j) {
163  V1(num_samples_obtained, j) = c[j];
164  }
165  proc_sampled_f_row[f_bv_max_global.proc].insert(f_bv_max_global.row);
166  proc_f_row_to_tmp_fs_row[f_bv_max_global.proc][f_bv_max_global.row] =
167  num_samples_obtained;
168  num_samples_obtained++;
169  }
170 
171  if (num_samples_obtained == 0)
172  {
173  f_bv_max_local.row_val = -DBL_MAX;
174  f_bv_max_local.proc = myid;
175  for (int i = 0; i < num_rows; ++i) {
176  double f_bv_val = fabs(Vo->item(i, 0));
177  if (f_bv_val > f_bv_max_local.row_val) {
178  f_bv_max_local.row_val = f_bv_val;
179  f_bv_max_local.row = i;
180  }
181  }
182  MPI_Allreduce(&f_bv_max_local, &f_bv_max_global, 1,
183  MaxRowType, RowInfoOp, MPI_COMM_WORLD);
184  // Now get the first sampled row of the basis
185  if (f_bv_max_global.proc == myid) {
186  for (int j = 0; j < num_basis_vectors; ++j) {
187  c[j] = Vo->item(f_bv_max_global.row, j);
188  }
189  }
190  MPI_Bcast(c.data(), num_basis_vectors, MPI_DOUBLE,
191  f_bv_max_global.proc, MPI_COMM_WORLD);
192  // Now add the first sampled row of the basis to tmp_fs.
193  for (int j = 0; j < num_basis_vectors; ++j) {
194  V1(0, j) = c[j];
195  }
196  proc_sampled_f_row[f_bv_max_global.proc].insert(f_bv_max_global.row);
197  proc_f_row_to_tmp_fs_row[f_bv_max_global.proc][f_bv_max_global.row] = 0;
198  num_samples_obtained++;
199  }
200 
201  if (num_samples_obtained < num_samples)
202  {
203  Vector A(num_rows, f_basis.distributed());
204  Vector noM(num_rows, f_basis.distributed());
205 
206  Matrix A0(num_basis_vectors - 1, num_basis_vectors - 1, false);
207  Matrix V1_last_col(num_basis_vectors - 1, 1, false);
208  Matrix tt(num_rows, num_basis_vectors - 1, f_basis.distributed());
209  Matrix tt1(num_rows, num_basis_vectors - 1, f_basis.distributed());
210  Matrix g1(tt.numRows(), tt.numColumns(), f_basis.distributed());
211  Matrix GG(tt1.numRows(), tt1.numColumns(), f_basis.distributed());
212  Vector ls_res_first_row(num_basis_vectors - 1, false);
213  Vector nV(num_basis_vectors, false);
214 
215  int start_idx = 2;
216  if (total_num_init_samples > 1)
217  {
218  start_idx += (total_num_init_samples - 1);
219  }
220  for (int i = start_idx; i <= num_samples; i++)
221  {
222  if (i <= num_basis_vectors)
223  {
224  A0.setSize(i - 1, i - 1);
225  V1_last_col.setSize(i - 1, 1);
226  for (int j = 0; j < num_samples_obtained; j++)
227  {
228  for (int k = 0; k < i - 1; k++)
229  {
230  A0(j, k) = V1(j, k);
231  }
232  V1_last_col(j, 0) = V1(j, i - 1);
233  }
234 
235  std::unique_ptr<Matrix> atA0 = V1_last_col.transposeMult(A0);
236  tt.setSize(num_rows, i - 1);
237  tt1.setSize(num_rows, i - 1);
238 
239  double ata = 0.0;
240  for (int j = 0; j < V1_last_col.numRows(); j++)
241  {
242  for (int k = 0; k < V1_last_col.numColumns(); k++)
243  {
244  ata += (V1_last_col(j, k) * V1_last_col(j, k));
245  }
246  }
247 
248  std::unique_ptr<Matrix> lhs = A0.transposeMult(A0);
249  lhs->inverse();
250  Matrix* rhs = NULL;
251  if (myid == 0)
252  {
253  rhs = new Matrix(num_rows + atA0->numRows(), i - 1, f_basis.distributed());
254  for (int k = 0; k < rhs->numColumns(); k++)
255  {
256  rhs->item(0, k) = atA0->item(0, k);
257  }
258  for (int j = 1; j < rhs->numRows(); j++)
259  {
260  for (int k = 0; k < rhs->numColumns(); k++)
261  {
262  rhs->item(j, k) = Vo->item(j - 1, k);
263  }
264  }
265  }
266  else
267  {
268  rhs = new Matrix(num_rows, i - 1, f_basis.distributed());
269  for (int j = 0; j < rhs->numRows(); j++)
270  {
271  for (int k = 0; k < rhs->numColumns(); k++)
272  {
273  rhs->item(j, k) = Vo->item(j, k);
274  }
275  }
276  }
277 
278  std::unique_ptr<Matrix> ls_res = rhs->mult(*lhs);
279  delete rhs;
280 
281  Matrix* c_T = NULL;
282  if (myid == 0)
283  {
284  c_T = new Matrix(ls_res->getData() + ls_res->numColumns(),
285  ls_res->numRows() - 1, ls_res->numColumns(), f_basis.distributed(), true);
286  }
287  else
288  {
289  c_T = new Matrix(ls_res->getData(),
290  ls_res->numRows(), ls_res->numColumns(), f_basis.distributed(), true);
291  }
292  std::unique_ptr<Matrix> Vo_first_i_columns = Vo->getFirstNColumns(i - 1);
293 
294  Vector b(num_rows, f_basis.distributed());
295  for (int j = 0; j < Vo_first_i_columns->numRows(); j++)
296  {
297  double tmp = 1.0;
298  for (int k = 0; k < Vo_first_i_columns->numColumns(); k++)
299  {
300  tmp += (Vo_first_i_columns->item(j, k) * c_T->item(j, k));
301  }
302  b(j) = tmp;
303  }
304 
305  for (int j = 0; j < num_rows; j++)
306  {
307  for (int zz = 0; zz < i - 1; zz++)
308  {
309  tt(j, zz) = Vo->item(j, zz) * Vo->item(j, i - 1);
310  }
311  }
312 
313  g1.setSize(tt.numRows(), tt.numColumns());
314  for (int j = 0; j < g1.numRows(); j++)
315  {
316  for (int k = 0; k < g1.numColumns(); k++)
317  {
318  g1(j, k) = atA0->item(0, k) + tt(j, k);
319  }
320  }
321 
322  Vector g3(num_rows, f_basis.distributed());
323  for (int j = 0; j < c_T->numRows(); j++)
324  {
325  double tmp = 0.0;
326  for (int k = 0; k < c_T->numColumns(); k++)
327  {
328  tmp += c_T->item(j, k) * g1(j, k);
329  }
330  g3(j) = tmp / b(j);
331  }
332 
333  for (int j = 0; j < num_rows; j++)
334  {
335  for (int zz = 0; zz < i - 1; zz++)
336  {
337  tt1(j, zz) = c_T->item(j, zz) * (Vo->item(j, i - 1) - g3(j));
338  }
339  }
340 
341  delete c_T;
342 
343  ls_res_first_row.setSize(ls_res->numColumns());
344  if (myid == 0) {
345  for (int j = 0; j < ls_res->numColumns(); ++j) {
346  c[j] = ls_res->item(0, j);
347  }
348  }
349  MPI_Bcast(c.data(), ls_res->numColumns(), MPI_DOUBLE,
350  0, MPI_COMM_WORLD);
351  for (int j = 0; j < ls_res->numColumns(); ++j) {
352  ls_res_first_row(j) = c[j];
353  }
354  GG.setSize(tt1.numRows(), tt1.numColumns());
355  for (int j = 0; j < GG.numRows(); j++)
356  {
357  for (int k = 0; k < GG.numColumns(); k++)
358  {
359  GG(j, k) = ls_res_first_row(k) + tt1(j, k);
360  }
361  }
362 
363  for (int j = 0; j < A.dim(); j++)
364  {
365  double tmp = 0.0;
366  for (int k = 0; k < g1.numColumns(); k++)
367  {
368  tmp += g1(j, k) * GG(j, k);
369  }
370  A(j) = std::max(0.0, ata +
371  (Vo->item(j, i - 1) * Vo->item(j, i - 1)) - tmp);
372  }
373 
374  nV.setSize(i);
375  for (int j = 0; j < i; j++)
376  {
377  nV(j) = 0.0;
378  for (int k = 0; k < num_samples_obtained; k++)
379  {
380  nV(j) += (V1(k, j) * V1(k, j));
381  }
382  }
383 
384  for (int j = 0; j < noM.dim(); j++)
385  {
386  noM(j) = 0.0;
387  for (int k = 0; k < i; k++)
388  {
389  noM(j) += std::log(nV(k) + (Vo->item(j, k) * Vo->item(j, k)));
390  }
391  }
392 
393  for (int j = 0; j < A.dim(); j++)
394  {
395  A(j) = std::log(fabs(A(j))) + std::log(b(j)) - noM(j);
396  }
397  }
398  else
399  {
400  Matrix curr_V1(V1.getData(), num_samples_obtained,
401  num_basis_vectors, false, true);
402  std::unique_ptr<Matrix> lhs = curr_V1.transposeMult(curr_V1);
403  lhs->inverse();
404 
405  std::unique_ptr<Matrix> ls_res = Vo->mult(*lhs);
406 
407  nV.setSize(num_basis_vectors);
408  for (int j = 0; j < num_basis_vectors; j++)
409  {
410  nV(j) = 0.0;
411  for (int k = 0; k < num_samples_obtained; k++)
412  {
413  nV(j) += (V1(k, j) * V1(k, j));
414  }
415  }
416 
417  for (int j = 0; j < noM.dim(); j++)
418  {
419  noM(j) = 0.0;
420  for (int k = 0; k < num_basis_vectors; k++)
421  {
422  noM(j) += std::log(nV(k) + (Vo->item(j, k) * Vo->item(j, k)));
423  }
424  }
425 
426  for (int j = 0; j < A.dim(); j++)
427  {
428  double tmp = 0.0;
429  for (int k = 0; k < Vo->numColumns(); k++)
430  {
431  tmp += Vo->item(j, k) * ls_res->item(j, k);
432  }
433  A(j) = std::log(1 + tmp) - noM(j);
434  }
435  }
436 
437  f_bv_max_local.row_val = -DBL_MAX;
438  f_bv_max_local.proc = myid;
439  for (int j = 0; j < num_rows; ++j) {
440  if (proc_f_row_to_tmp_fs_row[f_bv_max_local.proc].find(j) ==
441  proc_f_row_to_tmp_fs_row[f_bv_max_local.proc].end())
442  {
443  double f_bv_val = A(j);
444  if (f_bv_val > f_bv_max_local.row_val) {
445  f_bv_max_local.row_val = f_bv_val;
446  f_bv_max_local.row = j;
447  }
448  }
449  }
450  MPI_Allreduce(&f_bv_max_local, &f_bv_max_global, 1,
451  MaxRowType, RowInfoOp, MPI_COMM_WORLD);
452  // Now get the first sampled row of the basis
453  if (f_bv_max_global.proc == myid) {
454  for (int j = 0; j < num_basis_vectors; ++j) {
455  c[j] = Vo->item(f_bv_max_global.row, j);
456  }
457  }
458  MPI_Bcast(c.data(), num_basis_vectors, MPI_DOUBLE,
459  f_bv_max_global.proc, MPI_COMM_WORLD);
460  // Now add the first sampled row of the basis to tmp_fs.
461  for (int j = 0; j < num_basis_vectors; ++j) {
462  V1(num_samples_obtained, j) = c[j];
463  }
464  proc_sampled_f_row[f_bv_max_global.proc].insert(f_bv_max_global.row);
465  proc_f_row_to_tmp_fs_row[f_bv_max_global.proc][f_bv_max_global.row] =
466  num_samples_obtained;
467  num_samples_obtained++;
468  }
469  }
470 
471  // Fill f_sampled_row, and f_sampled_rows_per_proc. Unscramble tmp_fs into
472  // f_basis_sampled_inv.
473  int idx = 0;
474  for (int i = 0; i < num_procs; ++i) {
475  set<int>& this_proc_sampled_f_row = proc_sampled_f_row[i];
476  map<int, int>& this_proc_f_row_to_tmp_fs_row =
477  proc_f_row_to_tmp_fs_row[i];
478  f_sampled_rows_per_proc[i] = this_proc_sampled_f_row.size();
479  for (set<int>::iterator j = this_proc_sampled_f_row.begin();
480  j != this_proc_sampled_f_row.end(); ++j) {
481  int this_f_row = *j;
482  f_sampled_row[idx] = this_f_row;
483  int tmp_fs_row = this_proc_f_row_to_tmp_fs_row[this_f_row];
484  for (int col = 0; col < num_f_basis_cols; ++col) {
485  f_basis_sampled_inv(idx, col) = V1(tmp_fs_row, col);
486  }
487  ++idx;
488  }
489  }
490 
491  CAROM_ASSERT(num_samples == idx);
492 
493  // Now invert f_basis_sampled_inv, storing its transpose.
494  f_basis_sampled_inv.transposePseudoinverse();
495 
496  // Free the MPI_Datatype and MPI_Op.
497  MPI_Type_free(&MaxRowType);
498  MPI_Op_free(&RowInfoOp);
499 }
500 
501 }
void setSize(int num_rows, int num_cols)
Sets the number of rows and columns of the matrix and reallocates storage if needed....
Definition: Matrix.h:149
double * getData() const
Get the matrix data as a pointer.
Definition: Matrix.h:844
void transposePseudoinverse()
Computes the transposePseudoinverse of this.
Definition: Matrix.cpp:605
bool distributed() const
Returns true if the Matrix is distributed.
Definition: Matrix.h:183
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
std::unique_ptr< Matrix > qr_factorize() const
Computes and returns the Q of the QR factorization of this.
Definition: Matrix.cpp:846
int numDistributedRows() const
Returns the number of rows of the Matrix across all processors.
Definition: Matrix.h:211
std::unique_ptr< Matrix > getFirstNColumns(int n) const
Get the first N columns of a matrix.
Definition: Matrix.cpp:260
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::unique_ptr< Matrix > mult(const Matrix &other) const
Multiplies this Matrix with other and returns the product.
Definition: Matrix.h:273
std::unique_ptr< Matrix > transposeMult(const Matrix &other) const
Multiplies the transpose of this Matrix with other and returns the product.
Definition: Matrix.h:491
void setSize(int dim)
Sets the length of the vector and reallocates storage if needed. All values are initialized to zero.
Definition: Vector.h:217
int dim() const
Returns the dimension of the Vector on this processor.
Definition: Vector.h:251
void S_OPT(const Matrix &f_basis, int num_f_basis_vectors_used, std::vector< int > &f_sampled_row, std::vector< int > &f_sampled_rows_per_proc, Matrix &f_basis_sampled_inv, const int myid, const int num_procs, const int num_samples_req, std::vector< int > *init_samples, bool qr_factorize)
Computes the S_OPT algorithm on the given basis.
Definition: S_OPT.cpp:30