libROM  v1.0
Data-driven physical simulation library
DMD.cpp
1 
11 // Description: Implementation of the DMD algorithm.
12 
13 #include "DMD.h"
14 
15 #include "linalg/Matrix.h"
16 #include "linalg/Vector.h"
17 #include "linalg/scalapack_wrapper.h"
18 #include "utils/Utilities.h"
19 #include "utils/CSVDatabase.h"
20 #include "utils/HDFDatabase.h"
21 #include "mpi.h"
22 
23 #include <cstring>
24 
25 /* Use C++11 built-in shared pointers if available; else fallback to Boost. */
26 #if __cplusplus >= 201103L
27 #include <memory>
28 #else
29 #include <boost/shared_ptr.hpp>
30 #endif
31 
32 /* Use automatically detected Fortran name-mangling scheme */
33 #define zgetrf CAROM_FC_GLOBAL(zgetrf, ZGETRF)
34 #define zgetri CAROM_FC_GLOBAL(zgetri, ZGETRI)
35 
36 extern "C" {
37  // LU decomposition of a general matrix.
38  void zgetrf(int*, int*, double*, int*, int*, int*);
39 
40  // Generate inverse of a matrix given its LU decomposition.
41  void zgetri(int*, double*, int*, int*, double*, int*, int*);
42 }
43 
44 namespace CAROM {
45 
46 DMD::DMD(int dim, bool alt_output_basis, std::shared_ptr<Vector> state_offset)
47 {
48  CAROM_VERIFY(dim > 0);
49 
50  // Get the rank of this process, and the number of processors.
51  int mpi_init;
52  MPI_Initialized(&mpi_init);
53  if (mpi_init == 0) {
54  MPI_Init(nullptr, nullptr);
55  }
56 
57  MPI_Comm_rank(MPI_COMM_WORLD, &d_rank);
58  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
59  d_dim = dim;
60  d_trained = false;
61  d_init_projected = false;
62  d_alt_output_basis = alt_output_basis;
63  setOffset(state_offset, 0);
64 }
65 
66 DMD::DMD(int dim, double dt, bool alt_output_basis,
67  std::shared_ptr<Vector> state_offset)
68 {
69  CAROM_VERIFY(dim > 0);
70  CAROM_VERIFY(dt > 0.0);
71 
72  // Get the rank of this process, and the number of processors.
73  int mpi_init;
74  MPI_Initialized(&mpi_init);
75  if (mpi_init == 0) {
76  MPI_Init(nullptr, nullptr);
77  }
78 
79  MPI_Comm_rank(MPI_COMM_WORLD, &d_rank);
80  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
81  d_dim = dim;
82  d_dt = dt;
83  d_trained = false;
84  d_init_projected = false;
85  d_alt_output_basis = alt_output_basis;
86  setOffset(state_offset, 0);
87 }
88 
89 DMD::DMD(std::string base_file_name)
90 {
91  // Get the rank of this process, and the number of processors.
92  int mpi_init;
93  MPI_Initialized(&mpi_init);
94  if (mpi_init == 0) {
95  MPI_Init(nullptr, nullptr);
96  }
97 
98  MPI_Comm_rank(MPI_COMM_WORLD, &d_rank);
99  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
100  d_trained = true;
101  d_init_projected = true;
102 
103  load(base_file_name);
104 }
105 
106 DMD::DMD(std::vector<std::complex<double>> & eigs,
107  std::shared_ptr<Matrix> & phi_real,
108  std::shared_ptr<Matrix> & phi_imaginary, int k,
109  double dt, double t_offset, std::shared_ptr<Vector> & state_offset)
110 {
111  // Get the rank of this process, and the number of processors.
112  int mpi_init;
113  MPI_Initialized(&mpi_init);
114  if (mpi_init == 0) {
115  MPI_Init(nullptr, nullptr);
116  }
117 
118  MPI_Comm_rank(MPI_COMM_WORLD, &d_rank);
119  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
120  d_trained = true;
121  d_init_projected = false;
122 
123  d_eigs = eigs;
124  d_phi_real = phi_real;
125  d_phi_imaginary = phi_imaginary;
126  d_k = k;
127  d_dt = dt;
128  d_t_offset = t_offset;
129  setOffset(state_offset, 0);
130 }
131 
132 void DMD::setOffset(std::shared_ptr<Vector> & offset_vector, int order)
133 {
134  if (order == 0)
135  {
136  d_state_offset = offset_vector;
137  }
138 }
139 
140 void DMD::takeSample(double* u_in, double t)
141 {
142  CAROM_VERIFY(u_in != 0);
143  CAROM_VERIFY(t >= 0.0);
144  Vector* sample = new Vector(u_in, d_dim, true);
145 
146  double orig_t = t;
147  if (d_snapshots.empty())
148  {
149  d_t_offset = t;
150  t = 0.0;
151  }
152  else
153  {
154  t -= d_t_offset;
155  }
156 
157  // Erase any snapshots taken at the same or later time
158  while (!d_sampled_times.empty() && d_sampled_times.back() >= t)
159  {
160  if (d_rank == 0) std::cout << "Removing existing snapshot at time: " <<
161  d_t_offset + d_sampled_times.back() << std::endl;
162  d_snapshots.pop_back();
163  d_sampled_times.pop_back();
164  }
165 
166  if (d_snapshots.empty())
167  {
168  d_t_offset = orig_t;
169  t = 0.0;
170  }
171  else
172  {
173  CAROM_VERIFY(d_sampled_times.back() < t);
174  }
175 
176  d_snapshots.push_back(std::shared_ptr<Vector>(sample));
177  d_sampled_times.push_back(t);
178 }
179 
180 void DMD::train(double energy_fraction, const Matrix* W0, double linearity_tol)
181 {
182  std::unique_ptr<const Matrix> f_snapshots = getSnapshotMatrix();
183  CAROM_VERIFY(f_snapshots->numColumns() > 1);
184  CAROM_VERIFY(energy_fraction > 0 && energy_fraction <= 1);
185  d_energy_fraction = energy_fraction;
186  constructDMD(*f_snapshots, d_rank, d_num_procs, W0, linearity_tol);
187 }
188 
189 void DMD::train(int k, const Matrix* W0, double linearity_tol)
190 {
191  std::unique_ptr<const Matrix> f_snapshots = getSnapshotMatrix();
192  CAROM_VERIFY(f_snapshots->numColumns() > 1);
193  CAROM_VERIFY(k > 0 && k <= f_snapshots->numColumns() - 1);
194  d_energy_fraction = -1.0;
195  d_k = k;
196  constructDMD(*f_snapshots, d_rank, d_num_procs, W0, linearity_tol);
197 }
198 
199 std::pair<Matrix*, Matrix*>
201 {
202  CAROM_VERIFY(snapshots.numColumns() > 1);
203 
204  // TODO: Making two copies of the snapshot matrix has a lot of overhead.
205  // We need to figure out a way to do submatrix multiplication and to
206  // reimplement this algorithm using one snapshot matrix.
207  Matrix* f_snapshots_in = new Matrix(snapshots.numRows(),
208  snapshots.numColumns() - 1, snapshots.distributed());
209  Matrix* f_snapshots_out = new Matrix(snapshots.numRows(),
210  snapshots.numColumns() - 1, snapshots.distributed());
211 
212  // Break up snapshots into snapshots_in and snapshots_out
213  // snapshots_in = all columns of snapshots except last
214  // snapshots_out = all columns of snapshots except first
215  for (int i = 0; i < snapshots.numRows(); i++)
216  {
217  for (int j = 0; j < snapshots.numColumns() - 1; j++)
218  {
219  f_snapshots_in->item(i, j) = snapshots(i, j);
220  f_snapshots_out->item(i, j) = snapshots(i, j + 1);
221  if (d_state_offset)
222  {
223  f_snapshots_in->item(i, j) -= d_state_offset->item(i);
224  f_snapshots_out->item(i, j) -= d_state_offset->item(i);
225  }
226  }
227  }
228 
229  return std::pair<Matrix*,Matrix*>(f_snapshots_in, f_snapshots_out);
230 }
231 
232 void
233 DMD::computePhi(DMDInternal dmd_internal_obj)
234 {
235  // Calculate phi
236  if (d_alt_output_basis)
237  {
238  std::unique_ptr<Matrix> f_snapshots_out_mult_d_basis_right =
239  dmd_internal_obj.snapshots_out->mult(*dmd_internal_obj.basis_right);
240  std::unique_ptr<Matrix> f_snapshots_out_mult_d_basis_right_mult_d_S_inv =
241  f_snapshots_out_mult_d_basis_right->mult(*dmd_internal_obj.S_inv);
242  d_phi_real = f_snapshots_out_mult_d_basis_right_mult_d_S_inv->mult(
243  *dmd_internal_obj.eigenpair->ev_real);
244  d_phi_imaginary = f_snapshots_out_mult_d_basis_right_mult_d_S_inv->mult(
245  *dmd_internal_obj.eigenpair->ev_imaginary);
246  }
247  else
248  {
249  d_phi_real = dmd_internal_obj.basis->mult(*dmd_internal_obj.eigenpair->ev_real);
250  d_phi_imaginary = dmd_internal_obj.basis->mult(
251  *dmd_internal_obj.eigenpair->ev_imaginary);
252  }
253 }
254 
255 void
256 DMD::constructDMD(const Matrix & f_snapshots,
257  int d_rank,
258  int d_num_procs,
259  const Matrix* W0,
260  double linearity_tol)
261 {
262  std::pair<Matrix*, Matrix*> f_snapshot_pair = computeDMDSnapshotPair(
263  f_snapshots);
264  Matrix* f_snapshots_in = f_snapshot_pair.first;
265  Matrix* f_snapshots_out = f_snapshot_pair.second;
266 
267  int *row_offset = new int[d_num_procs + 1];
268  row_offset[d_num_procs] = f_snapshots_in->numDistributedRows();
269  row_offset[d_rank] = f_snapshots_in->numRows();
270 
271  CAROM_VERIFY(MPI_Allgather(MPI_IN_PLACE,
272  1,
273  MPI_INT,
274  row_offset,
275  1,
276  MPI_INT,
277  MPI_COMM_WORLD) == MPI_SUCCESS);
278  for (int i = d_num_procs - 1; i >= 0; i--) {
279  row_offset[i] = row_offset[i + 1] - row_offset[i];
280  }
281 
282  CAROM_VERIFY(row_offset[0] == 0);
283 
284  int d_blocksize = row_offset[d_num_procs] / d_num_procs;
285  if (row_offset[d_num_procs] % d_num_procs != 0) d_blocksize += 1;
286 
287  SLPK_Matrix svd_input;
288 
289  // Calculate svd of snapshots_in
290  initialize_matrix(&svd_input, f_snapshots_in->numColumns(),
291  f_snapshots_in->numDistributedRows(),
292  1, d_num_procs, d_blocksize, d_blocksize);
293 
294  for (int rank = 0; rank < d_num_procs; ++rank)
295  {
296  scatter_block(&svd_input, 1, row_offset[rank] + 1,
297  f_snapshots_in->getData(),
298  f_snapshots_in->numColumns(),
299  row_offset[rank + 1] - row_offset[rank], rank);
300  }
301 
302  std::unique_ptr<SVDManager> d_factorizer(new SVDManager);
303 
304  // This block does the actual ScaLAPACK call to do the factorization.
305  svd_init(d_factorizer.get(), &svd_input);
306  d_factorizer->dov = 1;
307  factorize(d_factorizer.get());
308  free_matrix_data(&svd_input);
309 
310  // Compute how many basis vectors we will actually use.
311  d_num_singular_vectors = std::min(f_snapshots_in->numColumns(),
312  f_snapshots_in->numDistributedRows());
313  for (int i = 0; i < d_num_singular_vectors; i++)
314  {
315  d_sv.push_back(d_factorizer->S[i]);
316  }
317 
318  if (d_energy_fraction != -1.0)
319  {
321  if (d_energy_fraction < 1.0)
322  {
323  double total_energy = 0.0;
324  for (int i = 0; i < d_num_singular_vectors; i++)
325  {
326  total_energy += d_factorizer->S[i];
327  }
328  double current_energy = 0.0;
329  for (int i = 0; i < d_num_singular_vectors; i++)
330  {
331  current_energy += d_factorizer->S[i];
332  if (current_energy / total_energy >= d_energy_fraction)
333  {
334  d_k = i + 1;
335  break;
336  }
337  }
338  }
339  }
340 
341  if (d_rank == 0) std::cout << "Using " << d_k << " basis vectors out of " <<
342  d_num_singular_vectors << "." << std::endl;
343 
344  // Allocate the appropriate matrices and gather their elements.
345  d_basis.reset(new Matrix(f_snapshots.numRows(), d_k,
346  f_snapshots.distributed()));
347  Matrix* d_S_inv = new Matrix(d_k, d_k, false);
348  Matrix* d_basis_right = new Matrix(f_snapshots_in->numColumns(), d_k, false);
349 
350  for (int d_rank = 0; d_rank < d_num_procs; ++d_rank) {
351  // V is computed in the transposed order so no reordering necessary.
352  gather_block(&d_basis->item(0, 0), d_factorizer->V,
353  1, row_offset[static_cast<unsigned>(d_rank)]+1,
354  d_k, row_offset[static_cast<unsigned>(d_rank) + 1] -
355  row_offset[static_cast<unsigned>(d_rank)],
356  d_rank);
357 
358  // gather_transposed_block does the same as gather_block, but transposes
359  // it; here, it is used to go from column-major to row-major order.
360  gather_transposed_block(&d_basis_right->item(0, 0), d_factorizer->U, 1, 1,
361  f_snapshots_in->numColumns(), d_k, d_rank);
362  }
363  delete[] row_offset;
364 
365  // Get inverse of singular values by multiplying by reciprocal.
366  for (int i = 0; i < d_k; ++i)
367  {
368  d_S_inv->item(i, i) = 1 / d_factorizer->S[static_cast<unsigned>(i)];
369  }
370 
371  std::unique_ptr<Matrix> Q;
372 
373  // If W0 is not null, we need to ensure it is in the basis of W.
374  if (W0 != NULL)
375  {
376  CAROM_VERIFY(W0->numRows() == d_basis->numRows());
377  CAROM_VERIFY(linearity_tol >= 0.0);
378  std::vector<int> lin_independent_cols_W;
379 
380  // Copy W0 and orthogonalize.
381  Matrix* d_basis_init = new Matrix(f_snapshots.numRows(), W0->numColumns(),
382  true);
383  for (int i = 0; i < d_basis_init->numRows(); i++)
384  {
385  for (int j = 0; j < W0->numColumns(); j++)
386  {
387  d_basis_init->item(i, j) = W0->item(i, j);
388  }
389  }
390  d_basis_init->orthogonalize();
391 
392  Vector W_col(f_snapshots.numRows(), f_snapshots.distributed());
393  Vector l(W0->numColumns(), true);
394  Vector W0l(f_snapshots.numRows(), f_snapshots.distributed());
395  // Find which columns of d_basis are linearly independent from W0
396  for (int j = 0; j < d_basis->numColumns(); j++)
397  {
398  // l = W0' * u
399  for (int i = 0; i < f_snapshots.numRows(); i++)
400  {
401  W_col.item(i) = d_basis->item(i, j);
402  }
403  d_basis_init->transposeMult(W_col, l);
404 
405  // W0l = W0 * l
406  d_basis_init->mult(l, W0l);
407 
408  // Compute k = sqrt(u.u - 2.0*l.l + basisl.basisl) which is ||u -
409  // basisl||_{2}. This is the error in the projection of u into the
410  // reduced order space and subsequent lifting back to the full
411  // order space.
412  double k = W_col.inner_product(W_col) - 2.0*l.inner_product(l) +
413  W0l.inner_product(W0l);
414  if (k <= 0)
415  {
416  k = 0;
417  }
418  else
419  {
420  k = sqrt(k);
421  }
422 
423  // Use k to see if the vector addressed by u is linearly dependent
424  // on the left singular vectors.
425  if (k >= linearity_tol)
426  {
427  lin_independent_cols_W.push_back(j);
428  }
429  }
430  delete d_basis_init;
431 
432  // Add the linearly independent columns of W to W0. Call this new basis W_new.
433  Matrix* d_basis_new = new Matrix(f_snapshots.numRows(),
434  W0->numColumns() + lin_independent_cols_W.size(), true);
435  for (int i = 0; i < d_basis_new->numRows(); i++)
436  {
437  for (int j = 0; j < W0->numColumns(); j++)
438  {
439  d_basis_new->item(i, j) = W0->item(i, j);
440  }
441  for (int j = 0; j < lin_independent_cols_W.size(); j++)
442  {
443  d_basis_new->item(i, j + W0->numColumns()) = d_basis->item(i,
444  lin_independent_cols_W[j]);
445  }
446  }
447 
448  // Orthogonalize W_new.
449  d_basis_new->orthogonalize();
450 
451  // Calculate Q = W* x W_new;
452  Q = d_basis->transposeMult(*d_basis_new);
453 
454  d_basis.reset(d_basis_new);
455 
456  d_k = d_basis_new->numColumns();
457  if (d_rank == 0) std::cout << "After adding W0, now using " << d_k <<
458  " basis vectors." << std::endl;
459  }
460 
461  // Calculate A_tilde = U_transpose * f_snapshots_out * V * inv(S)
462  std::unique_ptr<Matrix> d_basis_mult_f_snapshots_out = d_basis->transposeMult(
463  *f_snapshots_out);
464  std::unique_ptr<Matrix> d_basis_mult_f_snapshots_out_mult_d_basis_right =
465  d_basis_mult_f_snapshots_out->mult(*d_basis_right);
466  if (Q == NULL)
467  {
468  d_A_tilde = d_basis_mult_f_snapshots_out_mult_d_basis_right->mult(*d_S_inv);
469  }
470  else
471  {
472  std::unique_ptr<Matrix>
473  d_basis_mult_f_snapshots_out_mult_d_basis_right_mult_d_S_inv =
474  d_basis_mult_f_snapshots_out_mult_d_basis_right->mult(*d_S_inv);
475  d_A_tilde = d_basis_mult_f_snapshots_out_mult_d_basis_right_mult_d_S_inv->mult(
476  *Q);
477  }
478 
479  // Calculate the right eigenvalues/eigenvectors of A_tilde
481  d_eigs = eigenpair.eigs;
482 
483  DMDInternal dmd_internal = {f_snapshots_in, f_snapshots_out, d_basis.get(), d_basis_right, d_S_inv, &eigenpair};
484  computePhi(dmd_internal);
485 
486  Vector init(f_snapshots_in->numRows(), true);
487  for (int i = 0; i < init.dim(); i++)
488  {
489  init(i) = f_snapshots_in->item(i, 0);
490  }
491 
492  // Calculate pinv(d_phi) * initial_condition.
494 
495  d_trained = true;
496 
497  delete d_basis_right;
498  delete d_S_inv;
499  delete f_snapshots_in;
500  delete f_snapshots_out;
501  delete eigenpair.ev_real;
502  delete eigenpair.ev_imaginary;
503 
504  free_matrix_data(d_factorizer->U);
505  free_matrix_data(d_factorizer->V);
506  free(d_factorizer->U);
507  free(d_factorizer->V);
508  free(d_factorizer->S);
509 
510  release_context(&svd_input);
511 }
512 
513 void
514 DMD::projectInitialCondition(const Vector & init, double t_offset)
515 {
517  std::unique_ptr<Matrix> d_phi_real_squared_2 = d_phi_imaginary->transposeMult(
518  *d_phi_imaginary);
519  *d_phi_real_squared_inverse += *d_phi_real_squared_2;
520 
522  std::unique_ptr<Matrix> d_phi_imaginary_squared_2 =
523  d_phi_imaginary->transposeMult(*d_phi_real);
524  *d_phi_imaginary_squared_inverse -= *d_phi_imaginary_squared_2;
525 
526  const int dprs_row = d_phi_real_squared_inverse->numRows();
527  const int dprs_col = d_phi_real_squared_inverse->numColumns();
528  double* inverse_input = new double[dprs_row * dprs_col * 2];
529  for (int i = 0; i < d_phi_real_squared_inverse->numRows(); i++)
530  {
531  int k = 0;
532  for (int j = 0; j < d_phi_real_squared_inverse->numColumns() * 2; j++)
533  {
534  if (j % 2 == 0)
535  {
536  inverse_input[d_phi_real_squared_inverse->numColumns() * 2 * i + j] =
537  d_phi_real_squared_inverse->item(i, k);
538  }
539  else
540  {
541  inverse_input[d_phi_imaginary_squared_inverse->numColumns() * 2 * i + j] =
543  k++;
544  }
545  }
546  }
547 
548  // Call lapack routines to do the inversion.
549  // Set up some stuff the lapack routines need.
550  int info;
551  int mtx_size = d_phi_real_squared_inverse->numColumns();
552  int lwork = mtx_size*mtx_size*std::max(10,d_num_procs);
553  int* ipiv = new int [mtx_size];
554  double* work = new double [lwork];
555 
556  // Now call lapack to do the inversion.
557  zgetrf(&mtx_size, &mtx_size, inverse_input, &mtx_size, ipiv, &info);
558  zgetri(&mtx_size, inverse_input, &mtx_size, ipiv, work, &lwork, &info);
559 
560  for (int i = 0; i < d_phi_real_squared_inverse->numRows(); i++)
561  {
562  int k = 0;
563  for (int j = 0; j < d_phi_real_squared_inverse->numColumns() * 2; j++)
564  {
565  if (j % 2 == 0)
566  {
567  d_phi_real_squared_inverse->item(i, k) =
568  inverse_input[d_phi_real_squared_inverse->numColumns() * 2 * i + j];
569  }
570  else
571  {
572  d_phi_imaginary_squared_inverse->item(i, k) =
573  inverse_input[d_phi_imaginary_squared_inverse->numColumns() * 2 * i + j];
574  k++;
575  }
576  }
577  }
578 
579  std::unique_ptr<Vector> rhs_real = d_phi_real->transposeMult(init);
580  std::unique_ptr<Vector> rhs_imaginary = d_phi_imaginary->transposeMult(init);
581 
582  std::unique_ptr<Vector> d_projected_init_real_1 =
583  d_phi_real_squared_inverse->mult(*rhs_real);
584  std::unique_ptr<Vector> d_projected_init_real_2 =
585  d_phi_imaginary_squared_inverse->mult(*rhs_imaginary);
586  d_projected_init_real = d_projected_init_real_1->plus(*d_projected_init_real_2);
587 
588  std::unique_ptr<Vector> d_projected_init_imaginary_1 =
589  d_phi_real_squared_inverse->mult(*rhs_imaginary);
590  std::unique_ptr<Vector> d_projected_init_imaginary_2 =
591  d_phi_imaginary_squared_inverse->mult(*rhs_real);
592  d_projected_init_imaginary = d_projected_init_imaginary_2->minus(
593  *d_projected_init_imaginary_1);
594 
595  delete [] inverse_input;
596  delete [] ipiv;
597  delete [] work;
598 
599  if (t_offset >= 0.0)
600  {
601  std::cout << "t_offset is updated from " << d_t_offset <<
602  " to " << t_offset << std::endl;
603  d_t_offset = t_offset;
604  }
605  d_init_projected = true;
606 }
607 
608 std::unique_ptr<Vector>
609 DMD::predict(double t, int deg)
610 {
611  CAROM_VERIFY(d_trained);
612  CAROM_VERIFY(d_init_projected);
613  CAROM_VERIFY(t >= 0.0);
614 
615  t -= d_t_offset;
616 
617  std::pair<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>> d_phi_pair =
618  phiMultEigs(t, deg);
619  std::shared_ptr<Matrix> d_phi_mult_eigs_real = d_phi_pair.first;
620  std::shared_ptr<Matrix> d_phi_mult_eigs_imaginary = d_phi_pair.second;
621 
622  std::unique_ptr<Vector> d_predicted_state_real_1 = d_phi_mult_eigs_real->mult(
624 
625  std::unique_ptr<Vector> d_predicted_state_real_2 =
626  d_phi_mult_eigs_imaginary->mult(*d_projected_init_imaginary);
627 
628  std::unique_ptr<Vector> d_predicted_state_real =
629  d_predicted_state_real_1->minus(*d_predicted_state_real_2);
630  addOffset(*d_predicted_state_real, t, deg);
631 
632  return d_predicted_state_real;
633 }
634 
635 void
636 DMD::addOffset(Vector & result, double t, int deg)
637 {
638  if (d_state_offset)
639  {
640  result += *d_state_offset;
641  }
642 }
643 
644 std::complex<double>
645 DMD::computeEigExp(std::complex<double> eig, double t)
646 {
647  return std::pow(eig, t / d_dt);
648 }
649 
650 std::pair<std::shared_ptr<Matrix>,std::shared_ptr<Matrix>>
651  DMD::phiMultEigs(double t, int deg)
652 {
653  Matrix* d_eigs_exp_real = new Matrix(d_k, d_k, false);
654  Matrix* d_eigs_exp_imaginary = new Matrix(d_k, d_k, false);
655 
656  for (int i = 0; i < d_k; i++)
657  {
658  std::complex<double> eig_exp = computeEigExp(d_eigs[i], t);
659  for (int k = 0; k < deg; ++k)
660  {
661  eig_exp *= d_eigs[i];
662  }
663  d_eigs_exp_real->item(i, i) = std::real(eig_exp);
664  d_eigs_exp_imaginary->item(i, i) = std::imag(eig_exp);
665  }
666 
667  std::shared_ptr<Matrix> d_phi_mult_eigs_real = d_phi_real->mult(
668  *d_eigs_exp_real);
669  std::unique_ptr<Matrix> d_phi_mult_eigs_real_2 = d_phi_imaginary->mult(
670  *d_eigs_exp_imaginary);
671  *d_phi_mult_eigs_real -= *d_phi_mult_eigs_real_2;
672  std::shared_ptr<Matrix> d_phi_mult_eigs_imaginary = d_phi_real->mult(
673  *d_eigs_exp_imaginary);
674  std::unique_ptr<Matrix> d_phi_mult_eigs_imaginary_2 = d_phi_imaginary->mult(
675  *d_eigs_exp_real);
676  *d_phi_mult_eigs_imaginary += *d_phi_mult_eigs_imaginary_2;
677 
678  delete d_eigs_exp_real;
679  delete d_eigs_exp_imaginary;
680 
681  return std::pair<std::shared_ptr<Matrix>,std::shared_ptr<Matrix>>
682  (d_phi_mult_eigs_real,
683  d_phi_mult_eigs_imaginary);
684 }
685 
686 double
688 {
689  return d_t_offset;
690 }
691 
692 std::unique_ptr<const Matrix>
694 {
696 }
697 
698 std::unique_ptr<const Matrix>
699 DMD::createSnapshotMatrix(const std::vector<std::shared_ptr<Vector>> &
700  snapshots)
701 {
702  CAROM_VERIFY(snapshots.size() > 0);
703  CAROM_VERIFY(snapshots[0]->dim() > 0);
704  for (int i = 0 ; i < snapshots.size() - 1; i++)
705  {
706  CAROM_VERIFY(snapshots[i]->dim() == snapshots[i + 1]->dim());
707  CAROM_VERIFY(snapshots[i]->distributed() == snapshots[i + 1]->distributed());
708  }
709 
710  Matrix* snapshot_mat = new Matrix(snapshots[0]->dim(), snapshots.size(),
711  snapshots[0]->distributed());
712 
713  for (int i = 0; i < snapshots[0]->dim(); i++)
714  {
715  for (int j = 0; j < snapshots.size(); j++)
716  {
717  snapshot_mat->item(i, j) = snapshots[j]->item(i);
718  }
719  }
720 
721  return std::unique_ptr<const Matrix>(snapshot_mat);
722 }
723 
724 void
725 DMD::load(std::string base_file_name)
726 {
727  CAROM_ASSERT(!base_file_name.empty());
728 
729  char tmp[100];
730  std::string full_file_name = base_file_name;
731  HDFDatabase database;
732  database.open(full_file_name, "r");
733 
734  sprintf(tmp, "dt");
735  database.getDouble(tmp, d_dt);
736 
737  sprintf(tmp, "t_offset");
738  database.getDouble(tmp, d_t_offset);
739 
740  sprintf(tmp, "k");
741  database.getInteger(tmp, d_k);
742 
743  sprintf(tmp, "num_eigs");
744  int num_eigs;
745  database.getInteger(tmp, num_eigs);
746 
747  std::vector<double> eigs_real;
748  std::vector<double> eigs_imag;
749 
750  sprintf(tmp, "eigs_real");
751  eigs_real.resize(num_eigs);
752  database.getDoubleArray(tmp, &eigs_real[0], num_eigs);
753 
754  sprintf(tmp, "eigs_imag");
755  eigs_imag.resize(num_eigs);
756  database.getDoubleArray(tmp, &eigs_imag[0], num_eigs);
757 
758  for (int i = 0; i < num_eigs; i++)
759  {
760  d_eigs.push_back(std::complex<double>(eigs_real[i], eigs_imag[i]));
761  }
762  database.close();
763 
764  full_file_name = base_file_name + "_basis";
765  d_basis.reset(new Matrix());
766  d_basis->read(full_file_name);
767 
768  full_file_name = base_file_name + "_A_tilde";
769  d_A_tilde.reset(new Matrix());
770  d_A_tilde->read(full_file_name);
771 
772  full_file_name = base_file_name + "_phi_real";
773  d_phi_real.reset(new Matrix());
774  d_phi_real->read(full_file_name);
775 
776  full_file_name = base_file_name + "_phi_imaginary";
777  d_phi_imaginary.reset(new Matrix());
778  d_phi_imaginary->read(full_file_name);
779 
780  full_file_name = base_file_name + "_phi_real_squared_inverse";
781  d_phi_real_squared_inverse.reset(new Matrix());
782  d_phi_real_squared_inverse->read(full_file_name);
783 
784  full_file_name = base_file_name + "_phi_imaginary_squared_inverse";
786  d_phi_imaginary_squared_inverse->read(full_file_name);
787 
788  full_file_name = base_file_name + "_projected_init_real";
789  d_projected_init_real.reset(new Vector());
790  d_projected_init_real->read(full_file_name);
791 
792  full_file_name = base_file_name + "_projected_init_imaginary";
793  d_projected_init_imaginary.reset(new Vector());
794  d_projected_init_imaginary->read(full_file_name);
795 
796  full_file_name = base_file_name + "_state_offset";
797  if (Utilities::file_exist(full_file_name + ".000000"))
798  {
799  d_state_offset.reset(new Vector());
800  d_state_offset->read(full_file_name);
801  }
802 
803  d_init_projected = true;
804  d_trained = true;
805 
806  MPI_Barrier(MPI_COMM_WORLD);
807 }
808 
809 void
810 DMD::load(const char* base_file_name)
811 {
812  load(std::string(base_file_name));
813 }
814 
815 void
816 DMD::save(std::string base_file_name)
817 {
818  CAROM_ASSERT(!base_file_name.empty());
819  CAROM_VERIFY(d_trained);
820 
821  if (d_rank == 0)
822  {
823  char tmp[100];
824  std::string full_file_name = base_file_name;
825  HDFDatabase database;
826  database.create(full_file_name);
827 
828  sprintf(tmp, "dt");
829  database.putDouble(tmp, d_dt);
830 
831  sprintf(tmp, "t_offset");
832  database.putDouble(tmp, d_t_offset);
833 
834  sprintf(tmp, "k");
835  database.putInteger(tmp, d_k);
836 
837  sprintf(tmp, "num_eigs");
838  database.putInteger(tmp, d_eigs.size());
839 
840  std::vector<double> eigs_real;
841  std::vector<double> eigs_imag;
842 
843  for (int i = 0; i < d_eigs.size(); i++)
844  {
845  eigs_real.push_back(d_eigs[i].real());
846  eigs_imag.push_back(d_eigs[i].imag());
847  }
848 
849  sprintf(tmp, "eigs_real");
850  database.putDoubleArray(tmp, &eigs_real[0], eigs_real.size());
851 
852  sprintf(tmp, "eigs_imag");
853  database.putDoubleArray(tmp, &eigs_imag[0], eigs_imag.size());
854  database.close();
855  }
856 
857  std::string full_file_name;
858 
859  if (d_basis != NULL)
860  {
861  full_file_name = base_file_name + "_basis";
862  d_basis->write(full_file_name);
863  }
864 
865  if (d_A_tilde != NULL)
866  {
867  full_file_name = base_file_name + "_A_tilde";
868  d_A_tilde->write(full_file_name);
869  }
870 
871  full_file_name = base_file_name + "_phi_real";
872  d_phi_real->write(full_file_name);
873 
874  full_file_name = base_file_name + "_phi_imaginary";
875  d_phi_imaginary->write(full_file_name);
876 
877  full_file_name = base_file_name + "_phi_real_squared_inverse";
878  d_phi_real_squared_inverse->write(full_file_name);
879 
880  full_file_name = base_file_name + "_phi_imaginary_squared_inverse";
881  d_phi_imaginary_squared_inverse->write(full_file_name);
882 
883  full_file_name = base_file_name + "_projected_init_real";
884  d_projected_init_real->write(full_file_name);
885 
886  full_file_name = base_file_name + "_projected_init_imaginary";
887  d_projected_init_imaginary->write(full_file_name);
888 
889  if (d_state_offset != nullptr)
890  {
891  full_file_name = base_file_name + "_state_offset";
892  d_state_offset->write(full_file_name);
893  }
894 
895  MPI_Barrier(MPI_COMM_WORLD);
896 }
897 
898 void
899 DMD::save(const char* base_file_name)
900 {
901  save(std::string(base_file_name));
902 }
903 
904 void
905 DMD::summary(std::string base_file_name)
906 {
907  if (d_rank == 0)
908  {
909  CSVDatabase* csv_db(new CSVDatabase);
910 
911  csv_db->putDoubleVector(base_file_name + "_singular_value.csv", d_sv,
913  csv_db->putComplexVector(base_file_name + "_eigenvalue.csv", d_eigs,
914  d_eigs.size());
915 
916  delete csv_db;
917  }
918 }
919 
920 }
void putDoubleVector(const std::string &file_name, const std::vector< double > &data, int nelements, const bool distributed=false) override
Writes a vector of doubles associated with the supplied filename to the currently open CSV database f...
Definition: CSVDatabase.cpp:94
void putComplexVector(const std::string &file_name, const std::vector< std::complex< double >> &data, int nelements)
Writes a vector of complex doubles associated with the supplied filename.
void constructDMD(const Matrix &f_snapshots, int rank, int num_procs, const Matrix *W0, double linearity_tol)
Construct the DMD object.
Definition: DMD.cpp:256
std::vector< std::complex< double > > d_eigs
A vector holding the complex eigenvalues of the eigenmodes.
Definition: DMD.h:452
virtual void addOffset(Vector &result, double t=0.0, int deg=0)
Add the appropriate offset when predicting the solution.
Definition: DMD.cpp:636
double d_energy_fraction
The energy fraction used to obtain the DMD modes.
Definition: DMD.h:397
int d_dim
The total dimension of the sample vector.
Definition: DMD.h:347
std::unique_ptr< Matrix > d_phi_real_squared_inverse
The real part of d_phi_squared_inverse.
Definition: DMD.h:432
virtual void train(double energy_fraction, const Matrix *W0=NULL, double linearity_tol=0.0)
Train the DMD model with energy fraction criterion.
Definition: DMD.cpp:180
bool d_init_projected
Whether the initial condition has been projected.
Definition: DMD.h:382
std::shared_ptr< Vector > d_projected_init_imaginary
The imaginary part of the projected initial condition.
Definition: DMD.h:447
virtual void takeSample(double *u_in, double t)
Sample the new state, u_in. Any samples in d_snapshots taken at the same or later time will be erased...
Definition: DMD.cpp:140
virtual std::pair< Matrix *, Matrix * > computeDMDSnapshotPair(const Matrix &snapshots)
Returns a pair of pointers to the minus and plus snapshot matrices.
Definition: DMD.cpp:200
std::shared_ptr< Matrix > d_phi_imaginary
The imaginary part of d_phi.
Definition: DMD.h:427
virtual std::complex< double > computeEigExp(std::complex< double > eig, double t)
Compute the appropriate exponential function when predicting the solution.
Definition: DMD.cpp:645
void summary(std::string base_file_name)
Output the DMD record in CSV files.
Definition: DMD.cpp:905
std::pair< std::shared_ptr< Matrix >, std::shared_ptr< Matrix > > phiMultEigs(double t, int deg=0)
Internal function to multiply d_phi with the eigenvalues.
Definition: DMD.cpp:651
std::unique_ptr< const Matrix > getSnapshotMatrix()
Get the snapshot matrix contained within d_snapshots.
Definition: DMD.cpp:693
double d_dt
The time step size between samples.
Definition: DMD.h:352
virtual void save(std::string base_file_name)
Save the object state to a file.
Definition: DMD.cpp:816
std::vector< double > d_sv
std::vector holding the signular values.
Definition: DMD.h:392
double getTimeOffset() const
Get the time offset contained within d_t_offset.
Definition: DMD.cpp:687
void projectInitialCondition(const Vector &init, double t_offset=-1.0)
Project new initial condition using d_phi. Calculate pinv(phi) x init, or more precisely,...
Definition: DMD.cpp:514
std::unique_ptr< const Matrix > createSnapshotMatrix(const std::vector< std::shared_ptr< Vector >> &snapshots)
Get the snapshot matrix contained within d_snapshots.
Definition: DMD.cpp:699
int d_num_singular_vectors
The maximum number of singular vectors.
Definition: DMD.h:387
std::vector< double > d_sampled_times
The stored times of each sample.
Definition: DMD.h:367
double d_t_offset
The time offset of the first sample.
Definition: DMD.h:357
std::vector< std::shared_ptr< Vector > > d_snapshots
std::vector holding the snapshots.
Definition: DMD.h:362
int d_rank
The rank of the process this object belongs to.
Definition: DMD.h:337
DMD()
Unimplemented default constructor.
std::shared_ptr< Matrix > d_basis
The left singular vector basis.
Definition: DMD.h:407
std::unique_ptr< Vector > predict(double t, int deg=0)
Predict state given a time. Uses the projected initial condition of the training dataset (the first c...
Definition: DMD.cpp:609
int d_num_procs
The number of processors being run on.
Definition: DMD.h:342
std::shared_ptr< Vector > d_state_offset
State offset in snapshot.
Definition: DMD.h:372
bool d_alt_output_basis
Whether to use the alternative basis for output, i.e. phi = U^(+)*V*Omega^(-1)*X.
Definition: DMD.h:412
virtual void load(std::string base_file_name)
Load the object state from a file.
Definition: DMD.cpp:725
virtual void setOffset(std::shared_ptr< Vector > &offset_vector, int order)
Set the offset of a certain order.
Definition: DMD.cpp:132
std::shared_ptr< Matrix > d_A_tilde
A_tilde.
Definition: DMD.h:417
bool d_trained
Whether the DMD has been trained or not.
Definition: DMD.h:377
std::shared_ptr< Vector > d_projected_init_real
The real part of the projected initial condition.
Definition: DMD.h:442
virtual void computePhi(DMDInternal dmd_internal_obj)
Compute phi.
Definition: DMD.cpp:233
std::shared_ptr< Matrix > d_phi_real
The real part of d_phi.
Definition: DMD.h:422
std::unique_ptr< Matrix > d_phi_imaginary_squared_inverse
The imaginary part of d_phi_squared_inverse.
Definition: DMD.h:437
int d_k
The number of columns used after obtaining the SVD decomposition.
Definition: DMD.h:402
void getInteger(const std::string &key, int &data)
Reads an integer associated with the supplied key from the currently open database file.
Definition: Database.h:181
void putDouble(const std::string &key, double data)
Writes a double associated with the supplied key to currently open database file.
Definition: Database.h:128
void putInteger(const std::string &key, int data)
Writes an integer associated with the supplied key to currently open database file.
Definition: Database.h:94
void getDouble(const std::string &key, double &data)
Reads a double associated with the supplied key from the currently open database file.
Definition: Database.h:215
virtual void getDoubleArray(const std::string &key, double *data, int nelements, const bool distributed=false) override
Reads an array of doubles associated with the supplied key from the currently open HDF5 database file...
virtual bool create(const std::string &file_name, const MPI_Comm comm=MPI_COMM_NULL) override
Creates a new HDF5 database file with the supplied name.
Definition: HDFDatabase.cpp:45
virtual bool open(const std::string &file_name, const std::string &type, const MPI_Comm comm=MPI_COMM_NULL) override
Opens an existing HDF5 database file with the supplied name.
Definition: HDFDatabase.cpp:76
virtual void putDoubleArray(const std::string &key, const double *const data, int nelements, const bool distributed=false) override
Writes an array of doubles associated with the supplied key to the currently open HDF5 database file.
virtual bool close() override
Closes the currently open HDF5 database file.
double * getData() const
Get the matrix data as a pointer.
Definition: Matrix.h:844
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
int numDistributedRows() const
Returns the number of rows of the Matrix across all processors.
Definition: Matrix.h:211
void orthogonalize(bool double_pass=false, double zero_tol=1.0e-15)
Orthonormalizes the matrix.
Definition: Matrix.cpp:1631
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
const double & item(int i) const
Const Vector member access.
Definition: Vector.h:468
double inner_product(const Vector &other) const
Inner product, reference form.
Definition: Vector.cpp:196
int dim() const
Returns the dimension of the Vector on this processor.
Definition: Vector.h:251
ComplexEigenPair NonSymmetricRightEigenSolve(const Matrix &A)
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix.
Definition: Matrix.cpp:2031
Matrix * ev_real
The real component of the eigenvectors in matrix form.
Definition: Matrix.h:1126
std::vector< std::complex< double > > eigs
The corresponding complex eigenvalues.
Definition: Matrix.h:1136
Matrix * ev_imaginary
The imaginary component of the eigenvectors in matrix form.
Definition: Matrix.h:1131
ComplexEigenPair * eigenpair
The resultant DMD eigenvalues and eigenvectors.
Definition: DMD.h:64
Matrix * snapshots_out
The snapshot vectors of output.
Definition: DMD.h:44
Matrix * basis
The left singular vector basis.
Definition: DMD.h:49
Matrix * basis_right
The right singular vector basis.
Definition: DMD.h:54
Matrix * S_inv
The inverse of singular values.
Definition: DMD.h:59
static bool file_exist(const std::string &filename)
Verifies if a file exists.
Definition: Utilities.cpp:70