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);
66 
73  virtual
74  const Matrix*
76 
83  virtual
84  const Matrix*
86 
92  virtual
93  const Vector*
95 
101  virtual
102  const Matrix*
104 
105 protected:
113  virtual
114  void
116  double* u) = 0;
117 
128  virtual
129  bool
131  double* u, bool add_without_increase = false);
132 
136  virtual
137  void
139 
150  void
151  constructQ(
152  double*& Q,
153  const Vector* l,
154  double k);
155 
170  bool
171  svd(
172  double* A,
173  Matrix*& U,
174  Matrix*& S,
175  Matrix*& V);
176 
187  virtual
188  void
190  const Matrix* A,
191  const Matrix* W,
192  const Matrix* sigma) = 0;
193 
206  virtual
207  void
209  const Vector* j,
210  const Matrix* A,
211  const Matrix* W,
212  Matrix* sigma) = 0;
213 
219  int
221  {
222  return d_num_samples;
223  }
224 
234  double
236  const Matrix* m);
237 
243 
248 
253 
257  int d_size;
258 
262  int d_rank;
263 
267  std::vector<int> d_proc_dims;
268 
272  long int d_total_dim;
273 
281 
286 
292 
296  std::string d_state_file_name;
297 
301  static const int COMMUNICATE_U;
302 
303 private:
307  IncrementalSVD();
308 
313  const IncrementalSVD& other);
314 
319  operator = (
320  const IncrementalSVD& rhs);
321 };
322 
323 }
324 
325 #endif
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...
double checkOrthogonality(const Matrix *m)
Computes and returns the orthogonality of m.
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.
virtual void addNewSample(const Vector *j, const Matrix *A, const Matrix *W, Matrix *sigma)=0
Add a new, unique sample to the SVD.
static const int COMMUNICATE_U
MPI message tag.
virtual void computeBasis()=0
Computes the current basis vectors.
virtual const Vector * getSingularValues()
Returns the singular values for the current time interval.
int d_max_basis_dimension
The maximum basis dimension.
long int d_total_dim
The total dimension of the system.
std::string d_state_file_name
The name of file to which state is saved or restored from.
virtual const Matrix * getTemporalBasis()
Returns the temporal basis vectors for the current time interval as a Matrix.
void constructQ(double *&Q, const Vector *l, double k)
Construct the matrix Q whose SVD is needed.
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.
int d_size
The total number of processors.
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 bool takeSample(double *u_in, bool add_without_increase=false)
Sample new state, u_in, at the given time.
virtual ~IncrementalSVD()
Destructor.
virtual void buildInitialSVD(double *u)=0
Constructs the first SVD.
bool d_skip_linearly_dependent
Whether to skip linearly dependent samples.
virtual const Matrix * getSnapshotMatrix()
Returns the snapshot matrix for the current time interval.
virtual const Matrix * getSpatialBasis()
Returns the basis vectors for the current time interval as a Matrix.
virtual void addLinearlyDependentSample(const Matrix *A, const Matrix *W, const Matrix *sigma)=0
Add a linearly dependent sample to the SVD.
int d_rank
The rank of the processor owning this object.
Definition: SVD.h:29
int d_num_samples
Number of samples stored for the current time interval.
Definition: SVD.h:151