libROM  v1.0
Data-driven physical simulation library
BasisWriter.cpp
1 
11 // Description: A class that writes basis vectors to a file.
12 
13 #include "BasisWriter.h"
14 #include "utils/HDFDatabase.h"
15 #include "utils/HDFDatabaseMPIO.h"
16 #include "Matrix.h"
17 #include "Vector.h"
18 #include "BasisGenerator.h"
19 #include "utils/Utilities.h"
20 
21 #include "mpi.h"
22 
23 namespace CAROM {
24 
25 BasisWriter::BasisWriter(
26  BasisGenerator & basis_generator,
27  const std::string& base_file_name,
28  Database::formats db_format) :
29  d_basis_generator(basis_generator),
30  full_file_name(""),
31  snap_file_name(""),
32  db_format_(db_format),
33  d_database(NULL),
34  d_snap_database(NULL)
35 {
36  CAROM_ASSERT(!base_file_name.empty());
37 
38  int mpi_init;
39  MPI_Initialized(&mpi_init);
40  int rank;
41  if (mpi_init) {
42  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
43  }
44  else {
45  rank = 0;
46  }
47 
48  full_file_name = base_file_name;
49  snap_file_name = base_file_name + "_snapshot";
50 
51  // create and open snapshot/basis database
52  if (db_format_ == Database::formats::HDF5)
53  {
54  d_snap_database = new HDFDatabase();
55  d_database = new HDFDatabase();
56  }
57  else if (db_format_ == Database::formats::HDF5_MPIO)
58  {
59  d_snap_database = new HDFDatabaseMPIO();
60  d_database = new HDFDatabaseMPIO();
61  }
62  else
63  CAROM_ERROR("BasisWriter only supports HDF5/HDF5_MPIO data format!\n");
64 }
65 
67 {
68  delete d_database;
69  delete d_snap_database;
70 }
71 
72 void
73 BasisWriter::writeBasis(const std::string& kind)
74 {
75  CAROM_ASSERT(kind == "basis" || kind == "snapshot");
76 
77  char tmp[100];
78 
79  if (kind == "basis") {
80  d_database->create(full_file_name, MPI_COMM_WORLD);
81 
82  std::shared_ptr<const Matrix> basis = d_basis_generator.getSpatialBasis();
83  /* spatial basis is always distributed */
84  CAROM_VERIFY(basis->distributed());
85  int num_rows = basis->numRows();
86  int nrows_infile = num_rows;
87  if (db_format_ == Database::formats::HDF5_MPIO)
88  MPI_Allreduce(MPI_IN_PLACE, &nrows_infile, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
89  sprintf(tmp, "spatial_basis_num_rows");
90  d_database->putInteger(tmp, nrows_infile);
91  int num_cols = basis->numColumns();
92  sprintf(tmp, "spatial_basis_num_cols");
93  d_database->putInteger(tmp, num_cols);
94  sprintf(tmp, "spatial_basis");
95  d_database->putDoubleArray(tmp, &basis->item(0, 0), num_rows*num_cols, true);
96 
97  if (d_basis_generator.updateRightSV()) {
98  std::shared_ptr<const Matrix> tbasis = d_basis_generator.getTemporalBasis();
99  /* temporal basis is always not distributed */
100  CAROM_VERIFY(!tbasis->distributed());
101  num_rows = tbasis->numRows();
102  sprintf(tmp, "temporal_basis_num_rows");
103  d_database->putInteger(tmp, num_rows);
104  num_cols = tbasis->numColumns();
105  sprintf(tmp, "temporal_basis_num_cols");
106  d_database->putInteger(tmp, num_cols);
107  sprintf(tmp, "temporal_basis");
108  d_database->putDoubleArray(tmp, &tbasis->item(0, 0), num_rows*num_cols);
109  }
110 
111  std::shared_ptr<const Vector> sv = d_basis_generator.getSingularValues();
112  /* singular values are always not distributed */
113  CAROM_VERIFY(!sv->distributed());
114  int sv_dim = sv->dim();
115  sprintf(tmp, "singular_value_size");
116  d_database->putInteger(tmp, sv_dim);
117  sprintf(tmp, "singular_value");
118  d_database->putDoubleArray(tmp, &sv->item(0), sv_dim);
119 
120  d_database->close();
121  }
122 
123  if (kind == "snapshot") {
124  d_snap_database->create(snap_file_name, MPI_COMM_WORLD);
125 
126  std::shared_ptr<const Matrix> snapshots =
127  d_basis_generator.getSnapshotMatrix();
128  /* snapshot matrix is always distributed */
129  CAROM_VERIFY(snapshots->distributed());
130  int num_rows = snapshots->numRows(); // d_dim
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);
136  int num_cols = snapshots->numColumns(); // d_num_samples
137  sprintf(tmp, "snapshot_matrix_num_cols");
138  d_snap_database->putInteger(tmp, num_cols);
139  sprintf(tmp, "snapshot_matrix");
140  d_snap_database->putDoubleArray(tmp, &snapshots->item(0,0), num_rows*num_cols,
141  true);
142 
143  d_snap_database->close();
144  }
145 }
146 
147 }
std::shared_ptr< const Matrix > getSnapshotMatrix()
Returns the snapshot matrix for the current time interval.
bool updateRightSV()
Check whether right basis vectors will be updated.
std::shared_ptr< const Matrix > getSpatialBasis()
Returns the basis vectors for the current time interval as a Matrix.
std::shared_ptr< const Vector > getSingularValues()
Returns the singular values for the current time interval as a Vector.
std::shared_ptr< const Matrix > getTemporalBasis()
Returns the temporal basis vectors for the current time interval as a Matrix.
void writeBasis(const std::string &kind="basis")
Write basis or state vectors generated by d_basis_generator.
Definition: BasisWriter.cpp:73
~BasisWriter()
Destructor.
Definition: BasisWriter.cpp:66
virtual bool create(const std::string &file_name, const MPI_Comm comm=MPI_COMM_NULL)
Creates a new database file with the supplied name.
Definition: Database.cpp:29
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.
Definition: Database.h:311
void putInteger(const std::string &key, int data)
Writes an integer associated with the supplied key to currently open database file.
Definition: Database.h:94
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.