13 #include "BasisWriter.h"
14 #include "utils/HDFDatabase.h"
15 #include "utils/HDFDatabaseMPIO.h"
18 #include "BasisGenerator.h"
19 #include "utils/Utilities.h"
25 BasisWriter::BasisWriter(
27 const std::string& base_file_name,
29 d_basis_generator(basis_generator),
32 db_format_(db_format),
36 CAROM_ASSERT(basis_generator != 0);
37 CAROM_ASSERT(!base_file_name.empty());
40 MPI_Initialized(&mpi_init);
43 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
49 full_file_name = base_file_name;
50 snap_file_name = base_file_name +
"_snapshot";
53 if (db_format_ == Database::formats::HDF5)
58 else if (db_format_ == Database::formats::HDF5_MPIO)
64 CAROM_ERROR(
"BasisWriter only supports HDF5/HDF5_MPIO data format!\n");
70 delete d_snap_database;
76 CAROM_ASSERT(kind ==
"basis" || kind ==
"snapshot");
80 if (kind ==
"basis") {
81 d_database->
create(full_file_name, MPI_COMM_WORLD);
86 int num_rows = basis->
numRows();
87 int nrows_infile = num_rows;
88 if (db_format_ == Database::formats::HDF5_MPIO)
89 MPI_Allreduce(MPI_IN_PLACE, &nrows_infile, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
90 sprintf(tmp,
"spatial_basis_num_rows");
93 sprintf(tmp,
"spatial_basis_num_cols");
95 sprintf(tmp,
"spatial_basis");
103 sprintf(tmp,
"temporal_basis_num_rows");
106 sprintf(tmp,
"temporal_basis_num_cols");
108 sprintf(tmp,
"temporal_basis");
115 int sv_dim = sv->
dim();
116 sprintf(tmp,
"singular_value_size");
118 sprintf(tmp,
"singular_value");
124 if (kind ==
"snapshot") {
125 d_snap_database->
create(snap_file_name, MPI_COMM_WORLD);
130 int num_rows = snapshots->
numRows();
131 int nrows_infile = num_rows;
132 if (db_format_ == Database::formats::HDF5_MPIO)
133 MPI_Allreduce(MPI_IN_PLACE, &nrows_infile, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
134 sprintf(tmp,
"snapshot_matrix_num_rows");
135 d_snap_database->
putInteger(tmp, nrows_infile);
137 sprintf(tmp,
"snapshot_matrix_num_cols");
139 sprintf(tmp,
"snapshot_matrix");
143 d_snap_database->
close();
const Matrix * getSnapshotMatrix()
Returns the snapshot matrix for the current time interval.
bool updateRightSV()
Check whether right basis vectors will be updated.
const Matrix * getTemporalBasis()
Returns the temporal basis vectors for the current time interval as a Matrix.
const Matrix * getSpatialBasis()
Returns the basis vectors for the current time interval as a Matrix.
const Vector * getSingularValues()
Returns the singular values for the current time interval as a Vector.
void writeBasis(const std::string &kind="basis")
Write basis or state vectors generated by d_basis_generator.
~BasisWriter()
Destructor.
virtual bool create(const std::string &file_name, const MPI_Comm comm=MPI_COMM_NULL)
Creates a new database file with the supplied name.
virtual bool close()=0
Closes the currently open database file.
formats
Implemented database file formats. Add to this enum each time a new database format is implemented.
void putInteger(const std::string &key, int data)
Writes an integer associated with the supplied key to currently open database file.
virtual void putDoubleArray(const std::string &key, const double *const data, int nelements, const bool distributed=false)=0
Writes an array of doubles associated with the supplied key to the currently open database file.
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...
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.
bool distributed() const
Returns true if the Vector is distributed.
const double & item(int i) const
Const Vector member access.
int dim() const
Returns the dimension of the Vector on this processor.