libROM  v1.0
Data-driven physical simulation library
ParametricDMDc.h
1 
11 // Description: Computes the ParametricDMDc algorithm to obtain DMDc model interpolant
12 // at desired parameter point by interpolation of DMDc 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. We extend this algorithm to DMDc.
19 
20 #ifndef included_ParametricDMDc_h
21 #define included_ParametricDMDc_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 
55 template <class T>
56 void getParametricDMDc(std::unique_ptr<T>& parametric_dmdc,
57  const std::vector<Vector>& parameter_points,
58  std::vector<T*>& dmdcs,
59  std::vector<std::shared_ptr<Matrix>> & controls,
60  std::shared_ptr<Matrix> & controls_interpolated,
61  const Vector & desired_point,
62  std::string rbf = "G",
63  std::string interp_method = "LS",
64  double closest_rbf_val = 0.9,
65  bool reorthogonalize_W = false)
66 {
67  CAROM_VERIFY(parameter_points.size() == dmdcs.size());
68  CAROM_VERIFY(dmdcs.size() > 1);
69  for (int i = 0; i < dmdcs.size() - 1; i++)
70  {
71  CAROM_VERIFY(dmdcs[i]->d_dt == dmdcs[i + 1]->d_dt);
72  CAROM_VERIFY(dmdcs[i]->d_k == dmdcs[i + 1]->d_k);
73  }
74  CAROM_VERIFY(closest_rbf_val >= 0.0 && closest_rbf_val <= 1.0);
75 
76  int mpi_init, rank;
77  MPI_Initialized(&mpi_init);
78  if (mpi_init == 0) {
79  MPI_Init(nullptr, nullptr);
80  }
81 
82  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
83 
84  std::vector<std::shared_ptr<CAROM::Matrix>> bases;
85  std::vector<std::shared_ptr<CAROM::Matrix>> A_tildes;
86  std::vector<std::shared_ptr<CAROM::Matrix>> B_tildes;
87 
88  for (int i = 0; i < dmdcs.size(); i++)
89  {
90  bases.push_back(dmdcs[i]->d_basis);
91  A_tildes.push_back(dmdcs[i]->d_A_tilde);
92  B_tildes.push_back(dmdcs[i]->d_B_tilde);
93  }
94 
95  int ref_point = getClosestPoint(parameter_points, desired_point);
96  std::vector<std::shared_ptr<CAROM::Matrix>> rotation_matrices =
97  obtainRotationMatrices(parameter_points, bases, ref_point);
98 
99  CAROM::MatrixInterpolator basis_interpolator(parameter_points,
100  rotation_matrices, bases, ref_point, "B", rbf, interp_method, closest_rbf_val);
101  std::shared_ptr<CAROM::Matrix> W = basis_interpolator.interpolate(desired_point,
102  reorthogonalize_W);
103 
104  CAROM::MatrixInterpolator A_tilde_interpolator(parameter_points,
105  rotation_matrices, A_tildes, ref_point, "R", rbf, interp_method,
106  closest_rbf_val);
107  std::shared_ptr<CAROM::Matrix> A_tilde = A_tilde_interpolator.interpolate(
108  desired_point);
109 
110  CAROM::MatrixInterpolator B_tilde_interpolator(parameter_points,
111  rotation_matrices, B_tildes, ref_point, "R", rbf, interp_method,
112  closest_rbf_val);
113  std::shared_ptr<CAROM::Matrix> B_tilde = B_tilde_interpolator.interpolate(
114  desired_point);
115 
116  CAROM::MatrixInterpolator control_interpolator(parameter_points,
117  rotation_matrices, controls, ref_point, "R", rbf, interp_method,
118  closest_rbf_val);
119 
120  controls_interpolated = control_interpolator.interpolate(desired_point);
121 
122  // Calculate the right eigenvalues/eigenvectors of A_tilde
123  ComplexEigenPair eigenpair = NonSymmetricRightEigenSolve(*A_tilde);
124  std::vector<std::complex<double>> eigs = eigenpair.eigs;
125 
126  // Calculate phi (phi = W * eigenvectors)
127  std::shared_ptr<CAROM::Matrix> phi_real = W->mult(*eigenpair.ev_real);
128  std::shared_ptr<CAROM::Matrix> phi_imaginary = W->mult(*eigenpair.ev_imaginary);
129 
130  parametric_dmdc.reset(new T(eigs, phi_real, phi_imaginary, B_tilde,
131  dmdcs[0]->d_k,dmdcs[0]->d_dt,
132  dmdcs[0]->d_t_offset, dmdcs[0]->d_state_offset, W));
133 
134  delete eigenpair.ev_real;
135  delete eigenpair.ev_imaginary;
136 }
137 
160 template <class T>
161 void getParametricDMDc(std::unique_ptr<T>& parametric_dmdc,
162  const std::vector<Vector>& parameter_points,
163  std::vector<std::string>& dmdc_paths,
164  std::vector<std::shared_ptr<Matrix>> & controls,
165  std::shared_ptr<Matrix> & controls_interpolated,
166  const Vector & desired_point,
167  std::string rbf = "G",
168  std::string interp_method = "LS",
169  double closest_rbf_val = 0.9,
170  bool reorthogonalize_W = false)
171 {
172  std::vector<T*> dmdcs;
173  for (int i = 0; i < dmdc_paths.size(); i++)
174  {
175  T* dmdc = new T(dmdc_paths[i]);
176  dmdcs.push_back(dmdc);
177  }
178 
179  getParametricDMDc(parametric_dmdc, parameter_points, dmdcs, controls,
180  controls_interpolated,
181  desired_point,
182  rbf, interp_method, closest_rbf_val,
183  reorthogonalize_W);
184 
185  for (int i = 0; i < dmdcs.size(); i++)
186  {
187  delete dmdcs[i];
188  }
189 }
190 
191 }
192 
193 #endif
std::shared_ptr< Matrix > interpolate(const Vector &point, bool orthogonalize=false)
Obtain the interpolated reduced matrix of the unsampled parameter point.
ComplexEigenPair NonSymmetricRightEigenSolve(const Matrix &A)
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix.
Definition: Matrix.cpp:2031
std::vector< std::shared_ptr< Matrix > > obtainRotationMatrices(const std::vector< Vector > &parameter_points, const std::vector< std::shared_ptr< Matrix >> &bases, int ref_point)
Obtain the rotation matrices for all the parameter points using the basis of the reference point.
void getParametricDMDc(std::unique_ptr< T > &parametric_dmdc, const std::vector< Vector > &parameter_points, std::vector< T * > &dmdcs, std::vector< std::shared_ptr< Matrix >> &controls, std::shared_ptr< Matrix > &controls_interpolated, const Vector &desired_point, std::string rbf="G", std::string interp_method="LS", double closest_rbf_val=0.9, bool reorthogonalize_W=false)
Constructor.
int getClosestPoint(const std::vector< const Vector * > &points, const Vector &test_point)
Get closest point to a test point among a group of points.
Definition: Vector.cpp:586
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