libROM  v1.0
Data-driven physical simulation library
IncrementalSVDFastUpdate.cpp
1 
11 // Description: The concrete implementation of the incremental SVD algorithm
12 // using Matthew Brand's "fast update" method.
13 
14 #include "IncrementalSVDFastUpdate.h"
15 #include "utils/HDFDatabase.h"
16 
17 #include "mpi.h"
18 
19 #include <cmath>
20 #include <limits>
21 
22 namespace CAROM {
23 
24 IncrementalSVDFastUpdate::IncrementalSVDFastUpdate(
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 IncrementalSVDFastUpdate::~IncrementalSVDFastUpdate()
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 void
90 IncrementalSVDFastUpdate::buildInitialSVD(
91  double* u)
92 {
93  CAROM_VERIFY(u != 0);
94 
95  // Build d_S for this new time interval.
96  d_S = new Vector(1, false);
97  Vector u_vec(u, d_dim, true);
98  double norm_u = u_vec.norm();
99  d_S->item(0) = norm_u;
100 
101  // Build d_Up for this new time interval.
102  d_Up = new Matrix(1, 1, false);
103  d_Up->item(0, 0) = 1.0;
104 
105  // Build d_U for this new time interval.
106  d_U = new Matrix(d_dim, 1, true);
107  for (int i = 0; i < d_dim; ++i) {
108  d_U->item(i, 0) = u[i]/norm_u;
109  }
110 
111  // Build d_W for this new time interval.
112  if (d_update_right_SV) {
113  d_W = new Matrix(1, 1, false);
114  d_W->item(0, 0) = 1.0;
115  }
116 
117  // We now have the first sample for the new time interval.
118  d_num_samples = 1;
119  d_num_rows_of_W = 1;
120 
121  // Compute the basis vectors for this time interval.
122  computeBasis();
123 
124 }
125 
126 void
127 IncrementalSVDFastUpdate::computeBasis()
128 {
129  d_basis = d_U->mult(d_Up);
130  if(d_update_right_SV)
131  {
132  delete d_basis_right;
133  d_basis_right = new Matrix(*d_W);
134  }
135 
136  if(d_rank == 0) {
137  std::cout << "d_num_samples = " << d_num_samples << "\n";
138  std::cout << "d_num_rows_of_W = " << d_num_rows_of_W << "\n";
139  std::cout << "d_singular_value_tol = " << d_singular_value_tol << "\n";
140  std::cout << "smallest SV = " << d_S->item(d_num_samples-1) << "\n";
141  if (d_num_samples > 1) {
142  std::cout << "next smallest SV = " << d_S->item(d_num_samples-2) << "\n";
143  }
144  }
145  // remove the smallest singular value if it is smaller than d_singular_value_tol
146  if ( (d_singular_value_tol != 0.0) &&
147  (d_S->item(d_num_samples-1) < d_singular_value_tol) &&
148  (d_num_samples != 1) ) {
149 
150  if (d_rank == 0) std::cout << "removing a small singular value!\n";
151 
152  Matrix* d_basis_new = new Matrix(d_dim, d_num_samples-1,
153  d_basis->distributed());
154  for (int row = 0; row < d_dim; ++row) {
155  for (int col = 0; col < d_num_samples-1; ++col) {
156  d_basis_new->item(row, col) = d_basis->item(row,col);
157  }
158  }
159  delete d_basis;
160  d_basis = d_basis_new;
161 
162  if (d_update_right_SV)
163  {
164  Matrix* d_basis_right_new = new Matrix(d_num_rows_of_W, d_num_samples-1,
165  d_basis_right->distributed());
166  for (int row = 0; row < d_num_rows_of_W; ++row) {
167  for (int col = 0; col < d_num_samples-1; ++col) {
168  d_basis_right_new->item(row, col) = d_basis_right->item(row,col);
169  }
170  }
171  delete d_basis_right;
172  d_basis_right = d_basis_right_new;
173  }
174  --d_num_samples;
175  }
176 
177  // Reorthogonalize if necessary.
178  if (fabs(checkOrthogonality(d_basis)) >
179  std::numeric_limits<double>::epsilon()*static_cast<double>(d_num_samples)) {
180  d_basis->orthogonalize();
181  }
182  if(d_update_right_SV)
183  {
184  if (fabs(checkOrthogonality(d_basis_right)) >
185  std::numeric_limits<double>::epsilon()*d_num_samples) {
186  d_basis_right->orthogonalize();
187  }
188  }
189 
190 }
191 
192 void
193 IncrementalSVDFastUpdate::addLinearlyDependentSample(
194  const Matrix* A,
195  const Matrix* W,
196  const Matrix* sigma)
197 {
198  CAROM_VERIFY(A != 0);
199  CAROM_VERIFY(sigma != 0);
200 
201  // Chop a row and a column off of A to form Amod. Also form
202  // d_S by chopping a row and a column off of sigma.
203  Matrix Amod(d_num_samples, d_num_samples, false);
204  for (int row = 0; row < d_num_samples; ++row) {
205  for (int col = 0; col < d_num_samples; ++col) {
206  Amod.item(row, col) = A->item(row, col);
207  if (row == col)
208  {
209  d_S->item(col) = sigma->item(row, col);
210  }
211  }
212  }
213 
214  // Multiply d_Up and Amod and put result into d_Up.
215  Matrix* Up_times_Amod = d_Up->mult(Amod);
216  delete d_Up;
217  d_Up = Up_times_Amod;
218 
219  Matrix* new_d_W;
220  if (d_update_right_SV) {
221  // The new d_W is the product of the current d_W extended by another row
222  // and column and W. The only new value in the extended version of d_W
223  // that is non-zero is the new lower right value and it is 1. We will
224  // construct this product without explicitly forming the extended version of
225  // d_W.
226  new_d_W = new Matrix(d_num_rows_of_W+1, d_num_samples, false);
227  for (int row = 0; row < d_num_rows_of_W; ++row) {
228  for (int col = 0; col < d_num_samples; ++col) {
229  double new_d_W_entry = 0.0;
230  for (int entry = 0; entry < d_num_samples; ++entry) {
231  new_d_W_entry += d_W->item(row, entry)*W->item(entry, col);
232  }
233  new_d_W->item(row, col) = new_d_W_entry;
234  }
235  }
236  for (int col = 0; col < d_num_samples; ++col) {
237  new_d_W->item(d_num_rows_of_W, col) = W->item(d_num_samples, col);
238  }
239  delete d_W;
240  d_W = new_d_W;
241  ++d_num_rows_of_W;
242  }
243 
244 }
245 
246 void
247 IncrementalSVDFastUpdate::addNewSample(
248  const Vector* j,
249  const Matrix* A,
250  const Matrix* W,
251  Matrix* sigma)
252 {
253  CAROM_VERIFY(j != 0);
254  CAROM_VERIFY(A != 0);
255  CAROM_VERIFY(sigma != 0);
256 
257  // Add j as a new column of d_U.
258  Matrix* newU = new Matrix(d_dim, d_num_samples+1, true);
259  for (int row = 0; row < d_dim; ++row) {
260  for (int col = 0; col < d_num_samples; ++col) {
261  newU->item(row, col) = d_U->item(row, col);
262  }
263  newU->item(row, d_num_samples) = j->item(row);
264  }
265  delete d_U;
266  d_U = newU;
267 
268  Matrix* new_d_W;
269  if (d_update_right_SV) {
270  // The new d_W is the product of the current d_W extended by another row
271  // and column and W. The only new value in the extended version of d_W
272  // that is non-zero is the new lower right value and it is 1. We will
273  // construct this product without explicitly forming the extended version of
274  // d_W.
275  new_d_W = new Matrix(d_num_rows_of_W+1, d_num_samples+1, false);
276  for (int row = 0; row < d_num_rows_of_W; ++row) {
277  for (int col = 0; col < d_num_samples+1; ++col) {
278  double new_d_W_entry = 0.0;
279  for (int entry = 0; entry < d_num_samples; ++entry) {
280  new_d_W_entry += d_W->item(row, entry)*W->item(entry, col);
281  }
282  new_d_W->item(row, col) = new_d_W_entry;
283  }
284  }
285  for (int col = 0; col < d_num_samples+1; ++col) {
286  new_d_W->item(d_num_rows_of_W, col) = W->item(d_num_samples, col);
287  }
288  delete d_W;
289  d_W = new_d_W;
290  }
291 
292  // The new d_Up is the product of the current d_Up extended by another row
293  // and column and A. The only new value in the extended version of d_Up
294  // that is non-zero is the new lower right value and it is 1. We will
295  // construct this product without explicitly forming the extended version of
296  // d_Up.
297  Matrix* new_d_Up = new Matrix(d_num_samples+1, d_num_samples+1, false);
298  for (int row = 0; row < d_num_samples; ++row) {
299  for (int col = 0; col < d_num_samples+1; ++col) {
300  double new_d_Up_entry = 0.0;
301  for (int entry = 0; entry < d_num_samples; ++entry) {
302  new_d_Up_entry += d_Up->item(row, entry)*A->item(entry, col);
303  }
304  new_d_Up->item(row, col) = new_d_Up_entry;
305  }
306  }
307  for (int col = 0; col < d_num_samples+1; ++col) {
308  new_d_Up->item(d_num_samples, col) = A->item(d_num_samples, col);
309  }
310  delete d_Up;
311  d_Up = new_d_Up;
312 
313  // d_S = sigma.
314  delete d_S;
315  int num_dim = std::min(sigma->numRows(), sigma->numColumns());
316  d_S = new Vector(num_dim, false);
317  for (int i = 0; i < num_dim; i++) {
318  d_S->item(i) = sigma->item(i,i);
319  }
320 
321  // We now have another sample.
322  ++d_num_samples;
323  ++d_num_rows_of_W;
324 
325  // Reorthogonalize if necessary.
326  long int max_U_dim;
327  if (d_num_samples > d_total_dim) {
328  max_U_dim = d_num_samples;
329  }
330  else {
331  max_U_dim = d_total_dim;
332  }
333 
334 }
335 
336 }