libROM  v1.0
Data-driven physical simulation library
ParametricDMD.h
1 
11 // Description: Computes the ParametricDMD algorithm to obtain DMD model interpolant
12 // at desired parameter point by interpolation of DMD models at training parameter points.
13 // The implemented dynamic mode decomposition algorithm is derived from
14 // Tu et. al's paper "On Dynamic Mode Decomposition: Theory and
15 // Applications": https://arxiv.org/abs/1312.0041
16 // The interpolation algorithm was adapted from "Gradient-based
17 // Constrained Optimization Using a Database of Linear Reduced-Order Models"
18 // by Y. Choi et al.
19 
20 #ifndef included_ParametricDMD_h
21 #define included_ParametricDMD_h
22 
23 #include "manifold_interp/MatrixInterpolator.h"
24 #include "linalg/Matrix.h"
25 #include "linalg/Vector.h"
26 #include "mpi.h"
27 
28 #include <complex>
29 
30 namespace CAROM {
31 
51 template <class T>
52 void getParametricDMD(T*& parametric_dmd,
53  std::vector<Vector*>& parameter_points,
54  std::vector<T*>& dmds,
55  Vector* desired_point,
56  std::string rbf = "G",
57  std::string interp_method = "LS",
58  double closest_rbf_val = 0.9,
59  bool reorthogonalize_W = false)
60 {
61  CAROM_VERIFY(parameter_points.size() == dmds.size());
62  CAROM_VERIFY(dmds.size() > 1);
63  for (int i = 0; i < dmds.size() - 1; i++)
64  {
65  CAROM_VERIFY(dmds[i]->d_dt == dmds[i + 1]->d_dt);
66  CAROM_VERIFY(dmds[i]->d_k == dmds[i + 1]->d_k);
67  }
68  CAROM_VERIFY(closest_rbf_val >= 0.0 && closest_rbf_val <= 1.0);
69 
70  int mpi_init, rank;
71  MPI_Initialized(&mpi_init);
72  if (mpi_init == 0) {
73  MPI_Init(nullptr, nullptr);
74  }
75 
76  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
77 
78  std::vector<CAROM::Matrix*> bases;
79  std::vector<CAROM::Matrix*> A_tildes;
80  for (int i = 0; i < dmds.size(); i++)
81  {
82  bases.push_back(dmds[i]->d_basis);
83  A_tildes.push_back(dmds[i]->d_A_tilde);
84  }
85 
86  int ref_point = getClosestPoint(parameter_points, desired_point);
87  std::vector<CAROM::Matrix*> rotation_matrices = obtainRotationMatrices(
88  parameter_points,
89  bases, ref_point);
90 
91  CAROM::MatrixInterpolator basis_interpolator(parameter_points,
92  rotation_matrices, bases, ref_point, "B", rbf, interp_method, closest_rbf_val);
93  CAROM::Matrix* W = basis_interpolator.interpolate(desired_point,
94  reorthogonalize_W);
95 
96  CAROM::MatrixInterpolator A_tilde_interpolator(parameter_points,
97  rotation_matrices, A_tildes, ref_point, "R", rbf, interp_method,
98  closest_rbf_val);
99  CAROM::Matrix* A_tilde = A_tilde_interpolator.interpolate(desired_point);
100 
101  // Calculate the right eigenvalues/eigenvectors of A_tilde
102  ComplexEigenPair eigenpair = NonSymmetricRightEigenSolve(A_tilde);
103  std::vector<std::complex<double>> eigs = eigenpair.eigs;
104 
105  // Calculate phi (phi = W * eigenvectors)
106  CAROM::Matrix* phi_real = W->mult(eigenpair.ev_real);
107  CAROM::Matrix* phi_imaginary = W->mult(eigenpair.ev_imaginary);
108 
109  parametric_dmd = new T(eigs, phi_real, phi_imaginary, dmds[0]->d_k,
110  dmds[0]->d_dt, dmds[0]->d_t_offset,
111  dmds[0]->d_state_offset);
112 
113  delete W;
114  delete A_tilde;
115  delete eigenpair.ev_real;
116  delete eigenpair.ev_imaginary;
117 
118  for (auto m : rotation_matrices)
119  delete m;
120 }
121 
140 template <class T>
141 void getParametricDMD(T*& parametric_dmd,
142  std::vector<Vector*>& parameter_points,
143  std::vector<std::string>& dmd_paths,
144  Vector* desired_point,
145  std::string rbf = "G",
146  std::string interp_method = "LS",
147  double closest_rbf_val = 0.9,
148  bool reorthogonalize_W = false)
149 {
150  std::vector<T*> dmds;
151  for (int i = 0; i < dmd_paths.size(); i++)
152  {
153  T* dmd = new T(dmd_paths[i]);
154  dmds.push_back(dmd);
155  }
156 
157  getParametricDMD(parametric_dmd, parameter_points, dmds, desired_point,
158  rbf, interp_method, closest_rbf_val,
159  reorthogonalize_W);
160  for (int i = 0; i < dmds.size(); i++)
161  {
162  delete dmds[i];
163  }
164 }
165 
166 }
167 
168 #endif
Matrix * mult(const Matrix &other) const
Multiplies this Matrix with other and returns the product, reference version.
Definition: Matrix.h:283
Matrix * interpolate(Vector *point, bool orthogonalize=false)
Obtain the interpolated reduced matrix of the unsampled parameter point.
void getParametricDMD(T *&parametric_dmd, std::vector< Vector * > &parameter_points, std::vector< T * > &dmds, Vector *desired_point, std::string rbf="G", std::string interp_method="LS", double closest_rbf_val=0.9, bool reorthogonalize_W=false)
Constructor.
Definition: ParametricDMD.h:52
struct ComplexEigenPair NonSymmetricRightEigenSolve(Matrix *A)
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix.
Definition: Matrix.cpp:2203
int getClosestPoint(std::vector< Vector * > &points, Vector *test_point)
Get closest point to a test point among a group of points.
Definition: Vector.cpp:647
std::vector< Matrix * > obtainRotationMatrices(std::vector< Vector * > parameter_points, std::vector< Matrix * > bases, int ref_point)
Obtain the rotation matrices for all the parameter points using the basis of the reference point.
Matrix * ev_real
The real component of the eigenvectors in matrix form.
Definition: Matrix.h:1376
std::vector< std::complex< double > > eigs
The corresponding complex eigenvalues.
Definition: Matrix.h:1386
Matrix * ev_imaginary
The imaginary component of the eigenvectors in matrix form.
Definition: Matrix.h:1381