14 #include "StaticSVD.h"
17 #include "linalg/scalapack_wrapper.h"
25 #define dgesdd CAROM_FC_GLOBAL(dgesdd, DGESDD)
28 void dgesdd(
char*,
int*,
int*,
double*,
int*,
29 double*,
double*,
int*,
double*,
int*,
30 double*,
int*,
int*,
int*);
38 d_samples(new SLPK_Matrix), d_factorizer(new SVDManager),
39 d_basis_is_current(false),
40 d_max_basis_dimension(options.max_basis_dimension),
41 d_singular_value_tol(options.singular_value_tol),
42 d_preserve_snapshot(options.static_svd_preserve_snapshot)
46 MPI_Initialized(&mpi_init);
48 MPI_Init(
nullptr,
nullptr);
51 MPI_Comm_rank(MPI_COMM_WORLD, &
d_rank);
103 bool add_without_increase)
105 CAROM_VERIFY(u_in != 0);
106 CAROM_NULL_USE(add_without_increase);
112 if (u_vec.
norm() == 0.0) {
140 std::shared_ptr<const Matrix>
155 std::shared_ptr<const Matrix>
170 std::shared_ptr<const Vector>
181 CAROM_ASSERT(
d_S != 0);
187 std::shared_ptr<const Matrix>
191 CAROM_ERROR(
"StaticSVD: snapshot matrix is modified after computeSVD."
192 " To preserve the snapshots, set Options::static_svd_preserve_snapshot to be true!\n");
197 int nrows =
d_dims[
static_cast<unsigned>(rank)];
198 int firstrow =
d_istarts[
static_cast<unsigned>(rank)] + 1;
212 SLPK_Matrix *snapshot_matrix = NULL;
222 snapshot_matrix =
new SLPK_Matrix;
232 transpose_submatrix(snapshot_matrix, 1,
233 d_istarts[
static_cast<unsigned>(rank)]+1,
241 if (d_preserve_snapshot)
243 snapshot_matrix =
new SLPK_Matrix;
263 sigma_cutoff = std::numeric_limits<int>::max();
276 int ncolumns = hard_cutoff < sigma_cutoff ? hard_cutoff : sigma_cutoff;
278 CAROM_VERIFY(ncolumns >= 0);
283 CAROM_VERIFY(ncolumns >= 0);
284 unsigned nc =
static_cast<unsigned>(ncolumns);
285 memset(&
d_S->item(0), 0, nc*
sizeof(
double));
294 1,
d_istarts[
static_cast<unsigned>(rank)]+1,
295 ncolumns,
d_dims[
static_cast<unsigned>(rank)], rank);
307 d_istarts[
static_cast<unsigned>(rank)]+1,
308 1,
d_dims[
static_cast<unsigned>(rank)],
315 for (
int i = 0; i < ncolumns; ++i)
321 printf(
"Distribution of sampler's A and U:\n");
324 MPI_Barrier(MPI_COMM_WORLD);
327 printf(
"Distribution of sampler's V:\n");
330 MPI_Barrier(MPI_COMM_WORLD);
333 printf(
"Computed singular values: ");
334 for (
int i = 0; i < ncolumns; ++i)
341 if ((transpose) || (d_preserve_snapshot))
343 free_matrix_data(snapshot_matrix);
344 delete snapshot_matrix;
362 MPI_Allgather(&
d_dim, 1, MPI_INT,
d_dims.data(), 1, MPI_INT, MPI_COMM_WORLD);
366 for (
unsigned i = 0; i < static_cast<unsigned>(
d_num_procs); ++i) {
std::shared_ptr< Matrix > d_basis_right
The globalized right basis vectors for the current time interval.
std::shared_ptr< Matrix > d_basis
The globalized basis vectors for the current time interval.
bool isFirstSample() const
Returns true if the next sample will result in a new time interval.
const int d_max_num_samples
The maximum number of samples.
int d_num_samples
Number of samples stored for the current time interval.
bool d_debug_algorithm
Flag to indicate if results of algorithm should be printed for debugging purposes.
std::shared_ptr< Matrix > d_snapshots
The globalized snapshot vectors for the current time interval.
std::shared_ptr< Matrix > d_W
The matrix U which is large.
const int d_dim
Dimension of the system.
std::shared_ptr< Vector > d_S
The vector S which is small.
double d_singular_value_tol
The tolerance for singular values below which to drop vectors.
virtual bool takeSample(double *u_in, bool add_without_increase=false)
Collect the new sample, u_in at the supplied time.
std::shared_ptr< const Matrix > getSpatialBasis() override
Returns the basis vectors for the current time interval as a Matrix.
int d_rank
The rank of the process this object belongs to.
int d_nprow
The number of processor rows in the grid.
std::shared_ptr< const Matrix > getTemporalBasis() override
Returns the temporal basis vectors for the current time interval as a Matrix.
virtual void computeSVD()
Gathers samples from all other processors to form complete sample of system and computes the SVD.
int d_npcol
The number of processor columns in the grid.
void delete_factorizer()
Delete the factorizer from ScaLAPACK.
void get_global_info()
Get the system's total row dimension and where my rows sit in the matrix.
void delete_samples()
Delete the samples from ScaLAPACK.
void broadcast_sample(const double *u_in)
Broadcast the sample to all processors.
bool d_basis_is_current
Flag to indicate if the basis vectors for the current time interval are up to date.
std::vector< int > d_istarts
The starting row (0-based) of the matrix that each process owns. Stored to avoid an MPI operation to ...
int d_total_dim
The total dimension of the system (row dimension)
std::shared_ptr< const Vector > getSingularValues() override
Returns the singular values for the current time interval.
int d_max_basis_dimension
The max number of basis vectors to return.
std::unique_ptr< SLPK_Matrix > d_samples
Current samples of the system.
std::unique_ptr< SVDManager > d_factorizer
Factorization manager object used to compute the SVD.
int d_blocksize
The block size used internally for computing the SVD.
bool isBasisCurrent()
Checks if the basis vectors are computed from the current snapshot.
int d_num_procs
The number of processors being run on.
std::vector< int > d_dims
The number of rows that each process owns. Stored to avoid an MPI operation to get this operation eve...
std::shared_ptr< const Matrix > getSnapshotMatrix() override
Returns the snapshot matrix for the current time interval.
double norm() const
Form the norm of this.