libROM  v1.0
Data-driven physical simulation library
IncrementalSVDBrand.cpp
1 
11 // Description: The concrete implementation of the incremental SVD algorithm
12 // using Matthew Brand's "fast update" method.
13 
14 #include "IncrementalSVDBrand.h"
15 #include "utils/HDFDatabase.h"
16 
17 #include "mpi.h"
18 
19 #include <cmath>
20 #include <limits>
21 
22 namespace CAROM {
23 
24 IncrementalSVDBrand::IncrementalSVDBrand(
25  Options options,
26  const std::string& basis_file_name) :
27  IncrementalSVD(
28  options,
29  basis_file_name),
30  d_Up(0),
31  d_singular_value_tol(options.singular_value_tol)
32 {
33  CAROM_VERIFY(options.singular_value_tol >= 0);
34 
35  // If the state of the SVD is to be restored, do it now. The base class,
36  // IncrementalSVD, has already opened the database and restored the state
37  // common to all incremental algorithms. This particular class must also
38  // read the state of d_Up and then compute the basis. If the database could
39  // not be found then we can not restore the state.
40  if (options.restore_state && d_state_database) {
41  // Read d_Up.
42  int num_rows;
43  d_state_database->getInteger("Up_num_rows", num_rows);
44  int num_cols;
45  d_state_database->getInteger("Up_num_cols", num_cols);
46  d_Up = new Matrix(num_rows, num_cols, true);
47  d_state_database->getDoubleArray("Up",
48  &d_Up->item(0, 0),
49  num_rows*num_cols);
50 
51  // Close and delete the database.
52  d_state_database->close();
53  delete d_state_database;
54 
55  // Compute the basis.
56  computeBasis();
57  }
58 }
59 
60 IncrementalSVDBrand::~IncrementalSVDBrand()
61 {
62  // If the state of the SVD is to be saved, then create the database now.
63  // The IncrementalSVD base class destructor will save d_S and d_U. This
64  // derived class must save its specific state data, d_Up.
65  //
66  // If there are multiple time intervals then saving and restoring the state
67  // does not make sense as there is not one, all encompassing, basis.
68  if (d_save_state && (!isFirstSample())) {
69  // Create state database file.
70  d_state_database = new HDFDatabase();
71  d_state_database->create(d_state_file_name);
72 
73  // Save d_Up.
74  int num_rows = d_Up->numRows();
75  d_state_database->putInteger("Up_num_rows", num_rows);
76  int num_cols = d_Up->numColumns();
77  d_state_database->putInteger("Up_num_cols", num_cols);
78  d_state_database->putDoubleArray("Up",
79  &d_Up->item(0, 0),
80  num_rows*num_cols);
81  }
82 
83  // Delete data members.
84  if (d_Up) {
85  delete d_Up;
86  }
87 }
88 
89 const Matrix*
90 IncrementalSVDBrand::getSpatialBasis()
91 {
92  updateSpatialBasis(); // WARNING: this is costly
93 
94  CAROM_ASSERT(d_basis != 0);
95  return d_basis;
96 }
97 
98 const Matrix*
99 IncrementalSVDBrand::getTemporalBasis()
100 {
101  updateTemporalBasis();
102 
103  CAROM_ASSERT(d_basis_right != 0);
104  return d_basis_right;
105 }
106 
107 void
108 IncrementalSVDBrand::buildInitialSVD(
109  double* u)
110 {
111  CAROM_VERIFY(u != 0);
112 
113  // Build d_S for this new time interval.
114  d_S = new Vector(1, false);
115  Vector u_vec(u, d_dim, true);
116  double norm_u = u_vec.norm();
117  d_S->item(0) = norm_u;
118 
119  // Build d_Up for this new time interval.
120  d_Up = new Matrix(1, 1, false);
121  d_Up->item(0, 0) = 1.0;
122 
123  // Build d_U for this new time interval.
124  d_U = new Matrix(d_dim, 1, true);
125  for (int i = 0; i < d_dim; ++i) {
126  d_U->item(i, 0) = u[i]/norm_u;
127  }
128 
129  // Build d_W for this new time interval.
130  if (d_update_right_SV) {
131  d_W = new Matrix(1, 1, false);
132  d_W->item(0, 0) = 1.0;
133  }
134 
135  // We now have the first sample for the new time interval.
136  d_num_samples = 1;
137  d_num_rows_of_W = 1;
138 
139 }
140 
141 bool
142 IncrementalSVDBrand::buildIncrementalSVD(
143  double* u, bool add_without_increase)
144 {
145  CAROM_VERIFY(u != 0);
146 
147  // Compute the projection error
148  // (accurate down to the machine precision)
149  Vector u_vec(u, d_dim, true);
150  Vector e_proj(u, d_dim, true);
151  e_proj -= *(d_U->mult(d_U->transposeMult(e_proj))); // Gram-Schmidt
152  e_proj -= *(d_U->mult(d_U->transposeMult(e_proj))); // Re-orthogonalization
153 
154  double k = e_proj.inner_product(e_proj);
155  if (k <= 0) {
156  if(d_rank == 0) printf("linearly dependent sample!\n");
157  k = 0;
158  }
159  else {
160  k = sqrt(k);
161  }
162 
163  // Use k to see if the vector addressed by u is linearly dependent
164  // on the left singular vectors.
165  bool linearly_dependent_sample;
166  if ( k < d_linearity_tol ) {
167  if(d_rank == 0) {
168  std::cout << "linearly dependent sample! k = " << k << "\n";
169  std::cout << "d_linearity_tol = " << d_linearity_tol << "\n";
170  }
171  k = 0;
172  linearly_dependent_sample = true;
173  } else if ( d_num_samples >= d_max_basis_dimension || add_without_increase ) {
174  k = 0;
175  linearly_dependent_sample = true;
176  }
177  // Check to see if the "number of samples" (in IncrementalSVD and
178  // its subclasses, d_num_samples appears to be equal to the number
179  // of columns of the left singular vectors) is greater than or equal
180  // to the dimension of snapshot vectors. If so, then the vector
181  // addressed by the pointer u must be linearly dependent on the left
182  // singular vectors.
183  else if (d_num_samples >= d_total_dim) {
184  linearly_dependent_sample = true;
185  }
186  else {
187  linearly_dependent_sample = false;
188  }
189 
190  // Create Q.
191  double* Q;
192  Vector* U_mult_u = new Vector(d_U->transposeMult(u_vec)->getData(),
193  d_num_samples,
194  false);
195  Vector* l = d_Up->transposeMult(U_mult_u);
196  constructQ(Q, l, k);
197  delete U_mult_u, l;
198 
199  // Now get the singular value decomposition of Q.
200  Matrix* A;
201  Matrix* W;
202  Matrix* sigma;
203  bool result = svd(Q, A, sigma, W);
204 
205  // Done with Q.
206  delete [] Q;
207 
208  // If the svd was successful then add the sample. Otherwise clean up and
209  // return.
210  if (result) {
211 
212  // We need to add the sample if it is not linearly dependent or if it is
213  // linearly dependent and we are not skipping linearly dependent samples.
214  if ((linearly_dependent_sample && !d_skip_linearly_dependent) ) {
215  // This sample is linearly dependent and we are not skipping linearly
216  // dependent samples.
217  if(d_rank == 0) std::cout << "adding linearly dependent sample!\n";
218  addLinearlyDependentSample(A, W, sigma);
219  delete sigma;
220  }
221  else if (!linearly_dependent_sample) {
222  // This sample is not linearly dependent.
223 
224  // Compute j
225  Vector* j = new Vector(e_proj.getData(), d_dim, false);
226  for (int i = 0; i < d_dim; ++i) {
227  j->item(i) /= k;
228  }
229 
230  // addNewSample will assign sigma to d_S hence it should not be
231  // deleted upon return.
232  addNewSample(j, A, W, sigma);
233  delete j;
234  }
235  delete A;
236  delete W;
237  }
238  else {
239  delete A;
240  delete W;
241  delete sigma;
242  }
243  return result;
244 }
245 
246 void
247 IncrementalSVDBrand::updateSpatialBasis()
248 {
249  d_basis = d_U->mult(d_Up);
250 
251  // remove the smallest singular value if it is smaller than d_singular_value_tol
252  if ( (d_singular_value_tol != 0.0) &&
253  (d_S->item(d_num_samples-1) < d_singular_value_tol) &&
254  (d_num_samples != 1) ) {
255 
256  if (d_rank == 0) {
257  std::cout <<
258  "removing a spatial basis corresponding to the small singular value!\n";
259  }
260 
261  Matrix* d_basis_new = new Matrix(d_dim, d_num_samples-1,
262  d_basis->distributed());
263  for (int row = 0; row < d_dim; ++row) {
264  for (int col = 0; col < d_num_samples-1; ++col) {
265  d_basis_new->item(row, col) = d_basis->item(row,col);
266  }
267  }
268  delete d_basis;
269  d_basis = d_basis_new;
270  }
271 
272  // Reorthogonalize if necessary.
273  // (not likely to be called anymore but left for safety)
274  if (fabs(checkOrthogonality(d_basis)) >
275  std::numeric_limits<double>::epsilon()*static_cast<double>(d_num_samples)) {
276  d_basis->orthogonalize();
277  }
278 
279 }
280 
281 void
282 IncrementalSVDBrand::updateTemporalBasis()
283 {
284  delete d_basis_right;
285  d_basis_right = new Matrix(*d_W);
286 
287  // remove the smallest singular value if it is smaller than d_singular_value_tol
288  if ( (d_singular_value_tol != 0.0) &&
289  (d_S->item(d_num_samples-1) < d_singular_value_tol) &&
290  (d_num_samples != 1) ) {
291 
292  if (d_rank == 0) {
293  std::cout <<
294  "removing a temporal basis corresponding to the small singular value!\n";
295  }
296 
297  Matrix* d_basis_right_new = new Matrix(d_num_rows_of_W, d_num_samples-1,
298  d_basis_right->distributed());
299  for (int row = 0; row < d_num_rows_of_W; ++row) {
300  for (int col = 0; col < d_num_samples-1; ++col) {
301  d_basis_right_new->item(row, col) = d_basis_right->item(row,col);
302  }
303  }
304  delete d_basis_right;
305  d_basis_right = d_basis_right_new;
306  }
307 
308  // Reorthogonalize if necessary.
309  // (not likely to be called anymore but left for safety)
310  if (fabs(checkOrthogonality(d_basis_right)) >
311  std::numeric_limits<double>::epsilon()*d_num_samples) {
312  d_basis_right->orthogonalize();
313  }
314 
315 }
316 
317 void
318 IncrementalSVDBrand::computeBasis()
319 {
320  if(d_rank == 0) {
321  std::cout << "d_num_samples = " << d_num_samples << "\n";
322  std::cout << "d_num_rows_of_W = " << d_num_rows_of_W << "\n";
323  std::cout << "d_singular_value_tol = " << d_singular_value_tol << "\n";
324  std::cout << "smallest SV = " << d_S->item(d_num_samples-1) << "\n";
325  if (d_num_samples > 1) {
326  std::cout << "next smallest SV = " << d_S->item(d_num_samples-2) << "\n";
327  }
328  }
329 
330  updateSpatialBasis();
331  if (d_update_right_SV)
332  {
333  updateTemporalBasis();
334  }
335 
336  // remove the smallest singular value if it is smaller than d_singular_value_tol
337  if ( (d_singular_value_tol != 0.0) &&
338  (d_S->item(d_num_samples-1) < d_singular_value_tol) &&
339  (d_num_samples != 1) ) {
340 
341  --d_num_samples;
342  }
343 }
344 
345 void
346 IncrementalSVDBrand::addLinearlyDependentSample(
347  const Matrix* A,
348  const Matrix* W,
349  const Matrix* sigma)
350 {
351  CAROM_VERIFY(A != 0);
352  CAROM_VERIFY(sigma != 0);
353 
354  // Chop a row and a column off of A to form Amod. Also form
355  // d_S by chopping a row and a column off of sigma.
356  Matrix Amod(d_num_samples, d_num_samples, false);
357  for (int row = 0; row < d_num_samples; ++row) {
358  for (int col = 0; col < d_num_samples; ++col) {
359  Amod.item(row, col) = A->item(row, col);
360  if (row == col)
361  {
362  d_S->item(col) = sigma->item(row, col);
363  }
364  }
365  }
366 
367  // Multiply d_Up and Amod and put result into d_Up.
368  Matrix* Up_times_Amod = d_Up->mult(Amod);
369  delete d_Up;
370  d_Up = Up_times_Amod;
371 
372  Matrix* new_d_W;
373  if (d_update_right_SV) {
374  // The new d_W is the product of the current d_W extended by another row
375  // and column and W. The only new value in the extended version of d_W
376  // that is non-zero is the new lower right value and it is 1. We will
377  // construct this product without explicitly forming the extended version of
378  // d_W.
379  new_d_W = new Matrix(d_num_rows_of_W+1, d_num_samples, false);
380  for (int row = 0; row < d_num_rows_of_W; ++row) {
381  for (int col = 0; col < d_num_samples; ++col) {
382  double new_d_W_entry = 0.0;
383  for (int entry = 0; entry < d_num_samples; ++entry) {
384  new_d_W_entry += d_W->item(row, entry)*W->item(entry, col);
385  }
386  new_d_W->item(row, col) = new_d_W_entry;
387  }
388  }
389  for (int col = 0; col < d_num_samples; ++col) {
390  new_d_W->item(d_num_rows_of_W, col) = W->item(d_num_samples, col);
391  }
392  delete d_W;
393  d_W = new_d_W;
394  ++d_num_rows_of_W;
395  }
396 
397 }
398 
399 void
400 IncrementalSVDBrand::addNewSample(
401  const Vector* j,
402  const Matrix* A,
403  const Matrix* W,
404  Matrix* sigma)
405 {
406  CAROM_VERIFY(j != 0);
407  CAROM_VERIFY(A != 0);
408  CAROM_VERIFY(sigma != 0);
409 
410  // Add j as a new column of d_U.
411  Matrix* newU = new Matrix(d_dim, d_num_samples+1, true);
412  for (int row = 0; row < d_dim; ++row) {
413  for (int col = 0; col < d_num_samples; ++col) {
414  newU->item(row, col) = d_U->item(row, col);
415  }
416  newU->item(row, d_num_samples) = j->item(row);
417  }
418  delete d_U;
419  d_U = newU;
420 
421  Matrix* new_d_W;
422  if (d_update_right_SV) {
423  // The new d_W is the product of the current d_W extended by another row
424  // and column and W. The only new value in the extended version of d_W
425  // that is non-zero is the new lower right value and it is 1. We will
426  // construct this product without explicitly forming the extended version of
427  // d_W.
428  new_d_W = new Matrix(d_num_rows_of_W+1, d_num_samples+1, false);
429  for (int row = 0; row < d_num_rows_of_W; ++row) {
430  for (int col = 0; col < d_num_samples+1; ++col) {
431  double new_d_W_entry = 0.0;
432  for (int entry = 0; entry < d_num_samples; ++entry) {
433  new_d_W_entry += d_W->item(row, entry)*W->item(entry, col);
434  }
435  new_d_W->item(row, col) = new_d_W_entry;
436  }
437  }
438  for (int col = 0; col < d_num_samples+1; ++col) {
439  new_d_W->item(d_num_rows_of_W, col) = W->item(d_num_samples, col);
440  }
441  delete d_W;
442  d_W = new_d_W;
443  }
444 
445  // The new d_Up is the product of the current d_Up extended by another row
446  // and column and A. The only new value in the extended version of d_Up
447  // that is non-zero is the new lower right value and it is 1. We will
448  // construct this product without explicitly forming the extended version of
449  // d_Up.
450  Matrix* new_d_Up = new Matrix(d_num_samples+1, d_num_samples+1, false);
451  for (int row = 0; row < d_num_samples; ++row) {
452  for (int col = 0; col < d_num_samples+1; ++col) {
453  double new_d_Up_entry = 0.0;
454  for (int entry = 0; entry < d_num_samples; ++entry) {
455  new_d_Up_entry += d_Up->item(row, entry)*A->item(entry, col);
456  }
457  new_d_Up->item(row, col) = new_d_Up_entry;
458  }
459  }
460  for (int col = 0; col < d_num_samples+1; ++col) {
461  new_d_Up->item(d_num_samples, col) = A->item(d_num_samples, col);
462  }
463  delete d_Up;
464  d_Up = new_d_Up;
465 
466  // d_S = sigma.
467  delete d_S;
468  int num_dim = std::min(sigma->numRows(), sigma->numColumns());
469  d_S = new Vector(num_dim, false);
470  for (int i = 0; i < num_dim; i++) {
471  d_S->item(i) = sigma->item(i,i);
472  }
473 
474  // We now have another sample.
475  ++d_num_samples;
476  ++d_num_rows_of_W;
477 
478  // Reorthogonalize if necessary.
479  long int max_U_dim;
480  if (d_num_samples > d_total_dim) {
481  max_U_dim = d_num_samples;
482  }
483  else {
484  max_U_dim = d_total_dim;
485  }
486  if (fabs(checkOrthogonality(d_Up)) >
487  std::numeric_limits<double>::epsilon()*static_cast<double>(max_U_dim)) {
488  d_Up->orthogonalize();
489  }
490  if (fabs(checkOrthogonality(d_U)) >
491  std::numeric_limits<double>::epsilon()*static_cast<double>(max_U_dim)) {
492  d_U->orthogonalize(); // Will not be called, but just in case
493  }
494 
495  if(d_update_right_SV)
496  {
497  if (fabs(checkOrthogonality(d_W)) >
498  std::numeric_limits<double>::epsilon()*d_num_samples) {
499  d_W->orthogonalize();
500  }
501  }
502 
503 }
504 
505 }