libROM  v1.0
Data-driven physical simulation library
IncrementalSVD.h
1 
11 // Description: The abstract incremental SVD algorithm defines algorithm
12 // interface.
13 
14 #ifndef included_IncrementalSVD_h
15 #define included_IncrementalSVD_h
16 
17 #include "SVD.h"
18 #include "linalg/Options.h"
19 #include "utils/Database.h"
20 
21 namespace CAROM {
22 
27 class IncrementalSVD : public SVD
28 {
29 public:
42  Options options,
43  const std::string& basis_file_name);
44 
48  virtual
50 
61  virtual
62  bool
63  takeSample(
64  double* u_in,
65  bool add_without_increase = false) override;
66 
73  std::shared_ptr<const Matrix>
74  getSpatialBasis() override;
75 
82  std::shared_ptr<const Matrix>
83  getTemporalBasis() override;
84 
90  std::shared_ptr<const Vector>
91  getSingularValues() override;
92 
98  std::shared_ptr<const Matrix>
99  getSnapshotMatrix() override;
100 
101 protected:
109  virtual
110  void
112  double* u) = 0;
113 
124  virtual
125  bool
127  double* u, bool add_without_increase = false);
128 
132  virtual
133  void
135 
145  void
146  constructQ(
147  double*& Q,
148  const Vector & l,
149  double k);
150 
165  bool
166  svd(
167  double* A,
168  Matrix*& U,
169  Matrix*& S,
170  Matrix*& V);
171 
179  virtual
180  void
182  const Matrix & A,
183  const Matrix & W,
184  const Matrix & sigma) = 0;
185 
194  virtual
195  void
197  const Vector & j,
198  const Matrix & A,
199  const Matrix & W,
200  const Matrix & sigma) = 0;
201 
207  int
209  {
210  return d_num_samples;
211  }
212 
220  double
222  const Matrix & m);
223 
229 
234 
239 
243  int d_size;
244 
248  int d_rank;
249 
253  std::vector<int> d_proc_dims;
254 
258  long int d_total_dim;
259 
267 
272 
278 
282  std::string d_state_file_name;
283 
287  static const int COMMUNICATE_U;
288 
289 private:
293  IncrementalSVD();
294 
299  const IncrementalSVD& other);
300 
305  operator = (
306  const IncrementalSVD& rhs);
307 };
308 
309 }
310 
311 #endif
std::shared_ptr< const Matrix > getSpatialBasis() override
Returns the basis vectors for the current time interval as a Matrix.
bool d_save_state
If true the state of the SVD will be written to disk when the object is deleted. If there are multipl...
std::shared_ptr< const Matrix > getTemporalBasis() override
Returns the temporal basis vectors for the current time interval as a Matrix.
std::vector< int > d_proc_dims
The dimension of the system on each processor.
Database * d_state_database
Pointer to the database that will hold saved state data if the state is to be saved.
bool d_update_right_SV
Whether to update the right singular vectors.
double d_linearity_tol
Tolerance to determine whether or not a sample is linearly dependent.
static const int COMMUNICATE_U
MPI message tag.
virtual void computeBasis()=0
Computes the current basis vectors.
int d_max_basis_dimension
The maximum basis dimension.
long int d_total_dim
The total dimension of the system.
std::shared_ptr< const Vector > getSingularValues() override
Returns the singular values for the current time interval.
virtual void addLinearlyDependentSample(const Matrix &A, const Matrix &W, const Matrix &sigma)=0
Add a linearly dependent sample to the SVD.
std::string d_state_file_name
The name of file to which state is saved or restored from.
std::shared_ptr< const Matrix > getSnapshotMatrix() override
Returns the snapshot matrix for the current time interval.
virtual bool buildIncrementalSVD(double *u, bool add_without_increase=false)
Adds the new sampled the state vector, u, to the system.
int numSamples()
The number of samples stored.
virtual void addNewSample(const Vector &j, const Matrix &A, const Matrix &W, const Matrix &sigma)=0
Add a new, unique sample to the SVD.
int d_size
The total number of processors.
void constructQ(double *&Q, const Vector &l, double k)
Construct the matrix Q whose SVD is needed.
bool svd(double *A, Matrix *&U, Matrix *&S, Matrix *&V)
Given a matrix, A, returns 2 of the 3 components of its singular value decomposition....
virtual ~IncrementalSVD()
Destructor.
virtual void buildInitialSVD(double *u)=0
Constructs the first SVD.
bool d_skip_linearly_dependent
Whether to skip linearly dependent samples.
int d_rank
The rank of the processor owning this object.
double checkOrthogonality(const Matrix &m)
Computes and returns the orthogonality of m.
virtual bool takeSample(double *u_in, bool add_without_increase=false) override
Sample new state, u_in, at the given time.
Definition: SVD.h:29
int d_num_samples
Number of samples stored for the current time interval.
Definition: SVD.h:146