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"
26 #if __cplusplus >= 201103L
29 #include <boost/shared_ptr.hpp>
33 #define zgetrf CAROM_FC_GLOBAL(zgetrf, ZGETRF)
34 #define zgetri CAROM_FC_GLOBAL(zgetri, ZGETRI)
38 void zgetrf(
int*,
int*,
double*,
int*,
int*,
int*);
41 void zgetri(
int*,
double*,
int*,
int*,
double*,
int*,
int*);
46 DMD::DMD(
int dim,
bool alt_output_basis, std::shared_ptr<Vector> state_offset)
48 CAROM_VERIFY(dim > 0);
52 MPI_Initialized(&mpi_init);
54 MPI_Init(
nullptr,
nullptr);
57 MPI_Comm_rank(MPI_COMM_WORLD, &
d_rank);
66 DMD::DMD(
int dim,
double dt,
bool alt_output_basis,
67 std::shared_ptr<Vector> state_offset)
69 CAROM_VERIFY(dim > 0);
70 CAROM_VERIFY(dt > 0.0);
74 MPI_Initialized(&mpi_init);
76 MPI_Init(
nullptr,
nullptr);
79 MPI_Comm_rank(MPI_COMM_WORLD, &
d_rank);
93 MPI_Initialized(&mpi_init);
95 MPI_Init(
nullptr,
nullptr);
98 MPI_Comm_rank(MPI_COMM_WORLD, &
d_rank);
103 load(base_file_name);
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)
113 MPI_Initialized(&mpi_init);
115 MPI_Init(
nullptr,
nullptr);
118 MPI_Comm_rank(MPI_COMM_WORLD, &
d_rank);
142 CAROM_VERIFY(u_in != 0);
143 CAROM_VERIFY(t >= 0.0);
160 if (
d_rank == 0) std::cout <<
"Removing existing snapshot at time: " <<
176 d_snapshots.push_back(std::shared_ptr<Vector>(sample));
183 CAROM_VERIFY(f_snapshots->numColumns() > 1);
184 CAROM_VERIFY(energy_fraction > 0 && energy_fraction <= 1);
192 CAROM_VERIFY(f_snapshots->numColumns() > 1);
193 CAROM_VERIFY(k > 0 && k <= f_snapshots->numColumns() - 1);
199 std::pair<Matrix*, Matrix*>
215 for (
int i = 0; i < snapshots.
numRows(); i++)
217 for (
int j = 0; j < snapshots.
numColumns() - 1; j++)
219 f_snapshots_in->
item(i, j) = snapshots(i, j);
220 f_snapshots_out->
item(i, j) = snapshots(i, j + 1);
229 return std::pair<Matrix*,Matrix*>(f_snapshots_in, f_snapshots_out);
238 std::unique_ptr<Matrix> f_snapshots_out_mult_d_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(
244 d_phi_imaginary = f_snapshots_out_mult_d_basis_right_mult_d_S_inv->mult(
260 double linearity_tol)
264 Matrix* f_snapshots_in = f_snapshot_pair.first;
265 Matrix* f_snapshots_out = f_snapshot_pair.second;
271 CAROM_VERIFY(MPI_Allgather(MPI_IN_PLACE,
277 MPI_COMM_WORLD) == MPI_SUCCESS);
279 row_offset[i] = row_offset[i + 1] - row_offset[i];
282 CAROM_VERIFY(row_offset[0] == 0);
287 SLPK_Matrix svd_input;
290 initialize_matrix(&svd_input, f_snapshots_in->
numColumns(),
296 scatter_block(&svd_input, 1, row_offset[rank] + 1,
299 row_offset[rank + 1] - row_offset[rank], rank);
302 std::unique_ptr<SVDManager> d_factorizer(
new SVDManager);
305 svd_init(d_factorizer.get(), &svd_input);
306 d_factorizer->dov = 1;
307 factorize(d_factorizer.get());
308 free_matrix_data(&svd_input);
315 d_sv.push_back(d_factorizer->S[i]);
323 double total_energy = 0.0;
326 total_energy += d_factorizer->S[i];
328 double current_energy = 0.0;
331 current_energy += d_factorizer->S[i];
341 if (
d_rank == 0) std::cout <<
"Using " <<
d_k <<
" basis vectors out of " <<
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)],
360 gather_transposed_block(&d_basis_right->
item(0, 0), d_factorizer->U, 1, 1,
366 for (
int i = 0; i <
d_k; ++i)
368 d_S_inv->
item(i, i) = 1 / d_factorizer->S[
static_cast<unsigned>(i)];
371 std::unique_ptr<Matrix> Q;
377 CAROM_VERIFY(linearity_tol >= 0.0);
378 std::vector<int> lin_independent_cols_W;
383 for (
int i = 0; i < d_basis_init->
numRows(); i++)
387 d_basis_init->
item(i, j) = W0->
item(i, j);
396 for (
int j = 0; j <
d_basis->numColumns(); j++)
399 for (
int i = 0; i < f_snapshots.
numRows(); i++)
406 d_basis_init->
mult(l, W0l);
425 if (k >= linearity_tol)
427 lin_independent_cols_W.push_back(j);
434 W0->
numColumns() + lin_independent_cols_W.size(),
true);
435 for (
int i = 0; i < d_basis_new->
numRows(); i++)
439 d_basis_new->
item(i, j) = W0->
item(i, j);
441 for (
int j = 0; j < lin_independent_cols_W.size(); j++)
444 lin_independent_cols_W[j]);
452 Q =
d_basis->transposeMult(*d_basis_new);
457 if (
d_rank == 0) std::cout <<
"After adding W0, now using " <<
d_k <<
458 " basis vectors." << std::endl;
462 std::unique_ptr<Matrix> d_basis_mult_f_snapshots_out =
d_basis->transposeMult(
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);
468 d_A_tilde = d_basis_mult_f_snapshots_out_mult_d_basis_right->mult(*d_S_inv);
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(
483 DMDInternal dmd_internal = {f_snapshots_in, f_snapshots_out,
d_basis.get(), d_basis_right, d_S_inv, &eigenpair};
487 for (
int i = 0; i < init.
dim(); i++)
489 init(i) = f_snapshots_in->
item(i, 0);
497 delete d_basis_right;
499 delete f_snapshots_in;
500 delete f_snapshots_out;
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);
510 release_context(&svd_input);
517 std::unique_ptr<Matrix> d_phi_real_squared_2 =
d_phi_imaginary->transposeMult(
522 std::unique_ptr<Matrix> d_phi_imaginary_squared_2 =
528 double* inverse_input =
new double[dprs_row * dprs_col * 2];
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];
557 zgetrf(&mtx_size, &mtx_size, inverse_input, &mtx_size, ipiv, &info);
558 zgetri(&mtx_size, inverse_input, &mtx_size, ipiv, work, &lwork, &info);
579 std::unique_ptr<Vector> rhs_real =
d_phi_real->transposeMult(init);
580 std::unique_ptr<Vector> rhs_imaginary =
d_phi_imaginary->transposeMult(init);
582 std::unique_ptr<Vector> d_projected_init_real_1 =
584 std::unique_ptr<Vector> d_projected_init_real_2 =
588 std::unique_ptr<Vector> d_projected_init_imaginary_1 =
590 std::unique_ptr<Vector> d_projected_init_imaginary_2 =
593 *d_projected_init_imaginary_1);
595 delete [] inverse_input;
601 std::cout <<
"t_offset is updated from " <<
d_t_offset <<
602 " to " << t_offset << std::endl;
608 std::unique_ptr<Vector>
613 CAROM_VERIFY(t >= 0.0);
617 std::pair<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>> d_phi_pair =
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;
622 std::unique_ptr<Vector> d_predicted_state_real_1 = d_phi_mult_eigs_real->mult(
625 std::unique_ptr<Vector> d_predicted_state_real_2 =
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);
632 return d_predicted_state_real;
647 return std::pow(eig, t /
d_dt);
650 std::pair<std::shared_ptr<Matrix>,std::shared_ptr<Matrix>>
656 for (
int i = 0; i <
d_k; i++)
659 for (
int k = 0; k < deg; ++k)
663 d_eigs_exp_real->
item(i, i) = std::real(eig_exp);
664 d_eigs_exp_imaginary->
item(i, i) = std::imag(eig_exp);
667 std::shared_ptr<Matrix> d_phi_mult_eigs_real =
d_phi_real->mult(
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(
676 *d_phi_mult_eigs_imaginary += *d_phi_mult_eigs_imaginary_2;
678 delete d_eigs_exp_real;
679 delete d_eigs_exp_imaginary;
681 return std::pair<std::shared_ptr<Matrix>,std::shared_ptr<Matrix>>
682 (d_phi_mult_eigs_real,
683 d_phi_mult_eigs_imaginary);
692 std::unique_ptr<const Matrix>
698 std::unique_ptr<const Matrix>
702 CAROM_VERIFY(snapshots.size() > 0);
703 CAROM_VERIFY(snapshots[0]->dim() > 0);
704 for (
int i = 0 ; i < snapshots.size() - 1; i++)
706 CAROM_VERIFY(snapshots[i]->dim() == snapshots[i + 1]->dim());
707 CAROM_VERIFY(snapshots[i]->distributed() == snapshots[i + 1]->distributed());
710 Matrix* snapshot_mat =
new Matrix(snapshots[0]->dim(), snapshots.size(),
711 snapshots[0]->distributed());
713 for (
int i = 0; i < snapshots[0]->dim(); i++)
715 for (
int j = 0; j < snapshots.size(); j++)
717 snapshot_mat->
item(i, j) = snapshots[j]->item(i);
721 return std::unique_ptr<const Matrix>(snapshot_mat);
727 CAROM_ASSERT(!base_file_name.empty());
730 std::string full_file_name = base_file_name;
732 database.
open(full_file_name,
"r");
737 sprintf(tmp,
"t_offset");
743 sprintf(tmp,
"num_eigs");
747 std::vector<double> eigs_real;
748 std::vector<double> eigs_imag;
750 sprintf(tmp,
"eigs_real");
751 eigs_real.resize(num_eigs);
754 sprintf(tmp,
"eigs_imag");
755 eigs_imag.resize(num_eigs);
758 for (
int i = 0; i < num_eigs; i++)
760 d_eigs.push_back(std::complex<double>(eigs_real[i], eigs_imag[i]));
764 full_file_name = base_file_name +
"_basis";
768 full_file_name = base_file_name +
"_A_tilde";
772 full_file_name = base_file_name +
"_phi_real";
776 full_file_name = base_file_name +
"_phi_imaginary";
780 full_file_name = base_file_name +
"_phi_real_squared_inverse";
784 full_file_name = base_file_name +
"_phi_imaginary_squared_inverse";
788 full_file_name = base_file_name +
"_projected_init_real";
792 full_file_name = base_file_name +
"_projected_init_imaginary";
796 full_file_name = base_file_name +
"_state_offset";
806 MPI_Barrier(MPI_COMM_WORLD);
812 load(std::string(base_file_name));
818 CAROM_ASSERT(!base_file_name.empty());
824 std::string full_file_name = base_file_name;
826 database.
create(full_file_name);
831 sprintf(tmp,
"t_offset");
837 sprintf(tmp,
"num_eigs");
840 std::vector<double> eigs_real;
841 std::vector<double> eigs_imag;
843 for (
int i = 0; i <
d_eigs.size(); i++)
845 eigs_real.push_back(
d_eigs[i].real());
846 eigs_imag.push_back(
d_eigs[i].imag());
849 sprintf(tmp,
"eigs_real");
852 sprintf(tmp,
"eigs_imag");
857 std::string full_file_name;
861 full_file_name = base_file_name +
"_basis";
862 d_basis->write(full_file_name);
867 full_file_name = base_file_name +
"_A_tilde";
871 full_file_name = base_file_name +
"_phi_real";
874 full_file_name = base_file_name +
"_phi_imaginary";
877 full_file_name = base_file_name +
"_phi_real_squared_inverse";
880 full_file_name = base_file_name +
"_phi_imaginary_squared_inverse";
883 full_file_name = base_file_name +
"_projected_init_real";
886 full_file_name = base_file_name +
"_projected_init_imaginary";
891 full_file_name = base_file_name +
"_state_offset";
895 MPI_Barrier(MPI_COMM_WORLD);
901 save(std::string(base_file_name));
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...
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.
std::vector< std::complex< double > > d_eigs
A vector holding the complex eigenvalues of the eigenmodes.
virtual void addOffset(Vector &result, double t=0.0, int deg=0)
Add the appropriate offset when predicting the solution.
double d_energy_fraction
The energy fraction used to obtain the DMD modes.
int d_dim
The total dimension of the sample vector.
std::unique_ptr< Matrix > d_phi_real_squared_inverse
The real part of d_phi_squared_inverse.
virtual void train(double energy_fraction, const Matrix *W0=NULL, double linearity_tol=0.0)
Train the DMD model with energy fraction criterion.
bool d_init_projected
Whether the initial condition has been projected.
std::shared_ptr< Vector > d_projected_init_imaginary
The imaginary part of the projected initial condition.
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...
virtual std::pair< Matrix *, Matrix * > computeDMDSnapshotPair(const Matrix &snapshots)
Returns a pair of pointers to the minus and plus snapshot matrices.
std::shared_ptr< Matrix > d_phi_imaginary
The imaginary part of d_phi.
virtual std::complex< double > computeEigExp(std::complex< double > eig, double t)
Compute the appropriate exponential function when predicting the solution.
void summary(std::string base_file_name)
Output the DMD record in CSV files.
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.
std::unique_ptr< const Matrix > getSnapshotMatrix()
Get the snapshot matrix contained within d_snapshots.
double d_dt
The time step size between samples.
virtual void save(std::string base_file_name)
Save the object state to a file.
std::vector< double > d_sv
std::vector holding the signular values.
double getTimeOffset() const
Get the time offset contained within d_t_offset.
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,...
std::unique_ptr< const Matrix > createSnapshotMatrix(const std::vector< std::shared_ptr< Vector >> &snapshots)
Get the snapshot matrix contained within d_snapshots.
int d_num_singular_vectors
The maximum number of singular vectors.
std::vector< double > d_sampled_times
The stored times of each sample.
double d_t_offset
The time offset of the first sample.
std::vector< std::shared_ptr< Vector > > d_snapshots
std::vector holding the snapshots.
int d_rank
The rank of the process this object belongs to.
DMD()
Unimplemented default constructor.
std::shared_ptr< Matrix > d_basis
The left singular vector basis.
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...
int d_num_procs
The number of processors being run on.
std::shared_ptr< Vector > d_state_offset
State offset in snapshot.
bool d_alt_output_basis
Whether to use the alternative basis for output, i.e. phi = U^(+)*V*Omega^(-1)*X.
virtual void load(std::string base_file_name)
Load the object state from a file.
virtual void setOffset(std::shared_ptr< Vector > &offset_vector, int order)
Set the offset of a certain order.
std::shared_ptr< Matrix > d_A_tilde
A_tilde.
bool d_trained
Whether the DMD has been trained or not.
std::shared_ptr< Vector > d_projected_init_real
The real part of the projected initial condition.
virtual void computePhi(DMDInternal dmd_internal_obj)
Compute phi.
std::shared_ptr< Matrix > d_phi_real
The real part of d_phi.
std::unique_ptr< Matrix > d_phi_imaginary_squared_inverse
The imaginary part of d_phi_squared_inverse.
int d_k
The number of columns used after obtaining the SVD decomposition.
void getInteger(const std::string &key, int &data)
Reads an integer associated with the supplied key from the currently open database file.
void putDouble(const std::string &key, double data)
Writes a double associated with the supplied key to currently open database file.
void putInteger(const std::string &key, int data)
Writes an integer associated with the supplied key to currently open database file.
void getDouble(const std::string &key, double &data)
Reads a double associated with the supplied key from the currently open database file.
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.
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.
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.
bool distributed() const
Returns true if the Matrix is distributed.
int numColumns() const
Returns the number of columns in the Matrix. This method will return the same value from each process...
int numDistributedRows() const
Returns the number of rows of the Matrix across all processors.
void orthogonalize(bool double_pass=false, double zero_tol=1.0e-15)
Orthonormalizes the matrix.
const double & item(int row, int col) const
Const Matrix member access. Matrix data is stored in row-major format.
int numRows() const
Returns the number of rows of the Matrix on this processor.
std::unique_ptr< Matrix > mult(const Matrix &other) const
Multiplies this Matrix with other and returns the product.
std::unique_ptr< Matrix > transposeMult(const Matrix &other) const
Multiplies the transpose of this Matrix with other and returns the product.
const double & item(int i) const
Const Vector member access.
double inner_product(const Vector &other) const
Inner product, reference form.
int dim() const
Returns the dimension of the Vector on this processor.
ComplexEigenPair NonSymmetricRightEigenSolve(const Matrix &A)
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix.
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.
ComplexEigenPair * eigenpair
The resultant DMD eigenvalues and eigenvectors.
Matrix * snapshots_out
The snapshot vectors of output.
Matrix * basis
The left singular vector basis.
Matrix * basis_right
The right singular vector basis.
Matrix * S_inv
The inverse of singular values.
static bool file_exist(const std::string &filename)
Verifies if a file exists.