20 #ifndef included_ParametricDMDc_h
21 #define included_ParametricDMDc_h
23 #include "manifold_interp/MatrixInterpolator.h"
24 #include "linalg/Matrix.h"
25 #include "linalg/Vector.h"
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)
67 CAROM_VERIFY(parameter_points.size() == dmdcs.size());
68 CAROM_VERIFY(dmdcs.size() > 1);
69 for (
int i = 0; i < dmdcs.size() - 1; i++)
71 CAROM_VERIFY(dmdcs[i]->d_dt == dmdcs[i + 1]->d_dt);
72 CAROM_VERIFY(dmdcs[i]->d_k == dmdcs[i + 1]->d_k);
74 CAROM_VERIFY(closest_rbf_val >= 0.0 && closest_rbf_val <= 1.0);
77 MPI_Initialized(&mpi_init);
79 MPI_Init(
nullptr,
nullptr);
82 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
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;
88 for (
int i = 0; i < dmdcs.size(); i++)
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);
96 std::vector<std::shared_ptr<CAROM::Matrix>> rotation_matrices =
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,
105 rotation_matrices, A_tildes, ref_point,
"R", rbf, interp_method,
107 std::shared_ptr<CAROM::Matrix> A_tilde = A_tilde_interpolator.
interpolate(
111 rotation_matrices, B_tildes, ref_point,
"R", rbf, interp_method,
113 std::shared_ptr<CAROM::Matrix> B_tilde = B_tilde_interpolator.
interpolate(
117 rotation_matrices, controls, ref_point,
"R", rbf, interp_method,
120 controls_interpolated = control_interpolator.
interpolate(desired_point);
124 std::vector<std::complex<double>> eigs = eigenpair.
eigs;
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);
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));
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)
172 std::vector<T*> dmdcs;
173 for (
int i = 0; i < dmdc_paths.size(); i++)
175 T* dmdc =
new T(dmdc_paths[i]);
176 dmdcs.push_back(dmdc);
180 controls_interpolated,
182 rbf, interp_method, closest_rbf_val,
185 for (
int i = 0; i < dmdcs.size(); i++)
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.
std::vector< std::shared_ptr< Matrix > > obtainRotationMatrices(const std::vector< Vector > ¶meter_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 > ¶metric_dmdc, const std::vector< Vector > ¶meter_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.
Matrix * ev_real
The real component of the eigenvectors in matrix form.
std::vector< std::complex< double > > eigs
The corresponding complex eigenvalues.
Matrix * ev_imaginary
The imaginary component of the eigenvectors in matrix form.