libROM  v1.0
Data-driven physical simulation library
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
CAROM::DMDc Class Reference

#include <DMDc.h>

Public Member Functions

 DMDc (int dim, int dim_c, double dt, Vector *state_offset=NULL)
 Constructor. Basic DMDc with uniform time step size. More...
 
 DMDc (std::string base_file_name)
 Constructor. DMDc from saved models. More...
 
virtual ~DMDc ()
 Destroy the DMDc object.
 
virtual void setOffset (Vector *offset_vector)
 Set the state offset.
 
virtual void takeSample (double *u_in, double t, double *f_in, bool last_step=false)
 Sample the new state, u_in. Any samples in d_snapshots taken at the same or later time will be erased. More...
 
virtual void train (double energy_fraction, const Matrix *B=NULL)
 Train the DMDc model with energy fraction criterion. The control matrix B may be available and used in training. It is yet to be implemented and requires consideration in sparse matrix multiplication in general. (See Section III.B. in https://arxiv.org/pdf/1409.6358.pdf) In default, the control matrix is unknown, with B = NULL. More...
 
virtual void train (int k, const Matrix *B=NULL)
 Train the DMDc model with specified reduced dimension. The control matrix B may be available and used in training. It is yet to be implemented and requires consideration in sparse matrix multiplication in general. (See Section III.B. in https://arxiv.org/pdf/1409.6358.pdf) In default, the control matrix is unknown, with B = NULL. More...
 
void project (const Vector *init, const Matrix *controls, double t_offset=-1.0)
 Project U using d_phi, where U is the initial condition and the controls. Calculate pinv(phi) x U, or more precisely, (phi* x phi)^{-1} x phi* x U, where phi* is the conjugate transpose. More...
 
Vectorpredict (double t)
 Predict state given a time. Uses the projected initial condition of the training dataset (the first column). More...
 
double getTimeOffset () const
 Get the time offset contained within d_t_offset.
 
int getNumSamples () const
 Returns the number of samples taken. More...
 
int getDimension () const
 
const MatrixgetSnapshotMatrix ()
 Get the snapshot matrix contained within d_snapshots.
 
virtual void load (std::string base_file_name)
 Load the object state from a file. More...
 
void load (const char *base_file_name)
 Load the object state from a file. More...
 
virtual void save (std::string base_file_name)
 Save the object state to a file. More...
 
void save (const char *base_file_name)
 Save the object state to a file. More...
 
void summary (std::string base_file_name)
 Output the DMDc record in CSV files.
 

Protected Member Functions

 DMDc (int dim, int dim_c, Vector *state_offset=NULL)
 Constructor. Variant of DMDc with non-uniform time step size. More...
 
 DMDc (std::vector< std::complex< double >> eigs, Matrix *phi_real, Matrix *phi_imaginary, Matrix *B_tilde, int k, double dt, double t_offset, Vector *state_offset, Matrix *basis=nullptr)
 Constructor. Specified from DMDc components. More...
 
 DMDc ()
 Unimplemented default constructor.
 
 DMDc (const DMDc &other)
 Unimplemented copy constructor.
 
DMDcoperator= (const DMDc &rhs)
 Unimplemented assignment operator.
 
std::pair< Matrix *, Matrix * > phiMultEigs (double t)
 Internal function to multiply d_phi with the eigenvalues.
 
void constructDMDc (const Matrix *f_snapshots, const Matrix *f_controls, int rank, int num_procs, const Matrix *B)
 Construct the DMDc object.
 
virtual std::pair< Matrix *, Matrix * > computeDMDcSnapshotPair (const Matrix *snapshots, const Matrix *controls, const Matrix *B)
 Returns a pair of pointers to the minus and plus snapshot matrices.
 
virtual std::complex< double > computeEigExp (std::complex< double > eig, double t)
 Compute the appropriate exponential function when predicting the solution.
 
virtual void addOffset (Vector *&result)
 Add the state offset when predicting the solution.
 
const MatrixcreateSnapshotMatrix (std::vector< Vector * > snapshots)
 Get the snapshot matrix contained within d_snapshots.
 

Protected Attributes

int d_rank
 The rank of the process this object belongs to.
 
int d_num_procs
 The number of processors being run on.
 
int d_dim
 The total dimension of the sample vector.
 
int d_dim_c
 The total dimension of the control vector.
 
double d_dt = -1.0
 The time step size between samples.
 
double d_t_offset
 The time offset of the first sample.
 
std::vector< Vector * > d_snapshots
 std::vector holding the snapshots.
 
std::vector< Vector * > d_controls
 std::vector holding the controls.
 
std::vector< Vector * > d_sampled_times
 The stored times of each sample.
 
Vectord_state_offset = NULL
 State offset in snapshot.
 
bool d_trained
 Whether the DMDc has been trained or not.
 
bool d_init_projected
 Whether the initial condition has been projected.
 
int d_num_singular_vectors
 The maximum number of singular vectors.
 
std::vector< double > d_sv
 std::vector holding the signular values.
 
double d_energy_fraction
 The energy fraction used to obtain the DMDc modes.
 
int d_k
 The number of columns used after obtaining the SVD decomposition.
 
Matrixd_basis = NULL
 The left singular vector basis.
 
Matrixd_A_tilde = NULL
 A_tilde.
 
Matrixd_B_tilde = NULL
 B_tilde.
 
Matrixd_phi_real = NULL
 The real part of d_phi.
 
Matrixd_phi_imaginary = NULL
 The imaginary part of d_phi.
 
Matrixd_phi_real_squared_inverse = NULL
 The real part of d_phi_squared_inverse.
 
Matrixd_phi_imaginary_squared_inverse = NULL
 The imaginary part of d_phi_squared_inverse.
 
Vectord_projected_init_real = NULL
 The real part of the projected initial condition.
 
Vectord_projected_init_imaginary = NULL
 The imaginary part of the projected initial condition.
 
Matrixd_projected_controls_real = NULL
 The real part of the projected controls.
 
Matrixd_projected_controls_imaginary = NULL
 The imaginary part of the projected controls.
 
std::vector< std::complex< double > > d_eigs
 A vector holding the complex eigenvalues of the eigenmodes.
 

Friends

void getParametricDMDc (DMDc *&parametric_dmdc, std::vector< Vector * > &parameter_points, std::vector< DMDc * > &dmdcs, std::vector< Matrix * > controls, Matrix *&controls_interpolated, Vector *desired_point, std::string rbf, std::string interp_method, double closest_rbf_val, bool reorthogonalize_W)
 Obtain DMD model interpolant at desired parameter point by interpolation of DMD models from training parameter points. More...
 

Detailed Description

Class DMDc implements the DMDc algorithm on a given snapshot matrix.

Definition at line 32 of file DMDc.h.

Constructor & Destructor Documentation

◆ DMDc() [1/4]

CAROM::DMDc::DMDc ( int  dim,
int  dim_c,
double  dt,
Vector state_offset = NULL 
)

Constructor. Basic DMDc with uniform time step size.

Parameters
[in]dimThe full-order state dimension.
[in]dim_cThe control dimension.
[in]dtThe dt between samples.
[in]state_offsetThe state offset.

Definition at line 68 of file DMDc.cpp.

◆ DMDc() [2/4]

CAROM::DMDc::DMDc ( std::string  base_file_name)

Constructor. DMDc from saved models.

Parameters
[in]base_file_nameThe base part of the filename of the database to load when restarting from a save.

Definition at line 91 of file DMDc.cpp.

◆ DMDc() [3/4]

CAROM::DMDc::DMDc ( int  dim,
int  dim_c,
Vector state_offset = NULL 
)
protected

Constructor. Variant of DMDc with non-uniform time step size.

Parameters
[in]dimThe full-order state dimension.
[in]dim_cThe control dimension.
[in]state_offsetThe state offset.

Definition at line 47 of file DMDc.cpp.

◆ DMDc() [4/4]

CAROM::DMDc::DMDc ( std::vector< std::complex< double >>  eigs,
Matrix phi_real,
Matrix phi_imaginary,
Matrix B_tilde,
int  k,
double  dt,
double  t_offset,
Vector state_offset,
Matrix basis = nullptr 
)
protected

Constructor. Specified from DMDc components.

Parameters
[in]eigsd_eigs
[in]phi_reald_phi_real
[in]phi_imaginaryd_phi_imaginary
[in]B_tilded_B_tilde
[in]kd_k
[in]dtd_dt
[in]t_offsetd_t_offset
[in]state_offsetd_state_offset
[in]basisd_basis, set by DMDc class for offline stages. When interpolating a new DMDc, we enter the interpolated basis explicitly

Definition at line 108 of file DMDc.cpp.

Member Function Documentation

◆ getNumSamples()

int CAROM::DMDc::getNumSamples ( ) const
inline

Returns the number of samples taken.

Returns
The number of samples taken.

Definition at line 135 of file DMDc.h.

◆ load() [1/2]

void CAROM::DMDc::load ( const char *  base_file_name)

Load the object state from a file.

Parameters
[in]base_file_nameThe base part of the filename to load the database from.

Definition at line 945 of file DMDc.cpp.

◆ load() [2/2]

void CAROM::DMDc::load ( std::string  base_file_name)
virtual

Load the object state from a file.

Parameters
[in]base_file_nameThe base part of the filename to load the database from.

Definition at line 856 of file DMDc.cpp.

◆ predict()

Vector * CAROM::DMDc::predict ( double  t)

Predict state given a time. Uses the projected initial condition of the training dataset (the first column).

Parameters
[in]tThe time of the output state

Definition at line 717 of file DMDc.cpp.

◆ project()

void CAROM::DMDc::project ( const Vector init,
const Matrix controls,
double  t_offset = -1.0 
)

Project U using d_phi, where U is the initial condition and the controls. Calculate pinv(phi) x U, or more precisely, (phi* x phi)^{-1} x phi* x U, where phi* is the conjugate transpose.

Parameters
[in]initThe initial condition.
[in]controlsThe controls.
[in]t_offsetThe initial time offset.

Definition at line 589 of file DMDc.cpp.

◆ save() [1/2]

void CAROM::DMDc::save ( const char *  base_file_name)

Save the object state to a file.

Parameters
[in]base_file_nameThe base part of the filename to save the database to.

Definition at line 1040 of file DMDc.cpp.

◆ save() [2/2]

void CAROM::DMDc::save ( std::string  base_file_name)
virtual

Save the object state to a file.

Parameters
[in]base_file_nameThe base part of the filename to save the database to.

Definition at line 951 of file DMDc.cpp.

◆ takeSample()

void CAROM::DMDc::takeSample ( double *  u_in,
double  t,
double *  f_in,
bool  last_step = false 
)
virtual

Sample the new state, u_in. Any samples in d_snapshots taken at the same or later time will be erased.

Precondition
u_in != 0
t >= 0.0
Parameters
[in]u_inThe new state.
[in]tThe time of the newly sampled state.
[in]f_inThe control.
[in]last_stepWhether it is the last step.

Definition at line 162 of file DMDc.cpp.

◆ train() [1/2]

void CAROM::DMDc::train ( double  energy_fraction,
const Matrix B = NULL 
)
virtual

Train the DMDc model with energy fraction criterion. The control matrix B may be available and used in training. It is yet to be implemented and requires consideration in sparse matrix multiplication in general. (See Section III.B. in https://arxiv.org/pdf/1409.6358.pdf) In default, the control matrix is unknown, with B = NULL.

Parameters
[in]energy_fractionThe energy fraction to keep after doing SVD.
[in]B(Optional) The control matrix B.

Definition at line 212 of file DMDc.cpp.

◆ train() [2/2]

void CAROM::DMDc::train ( int  k,
const Matrix B = NULL 
)
virtual

Train the DMDc model with specified reduced dimension. The control matrix B may be available and used in training. It is yet to be implemented and requires consideration in sparse matrix multiplication in general. (See Section III.B. in https://arxiv.org/pdf/1409.6358.pdf) In default, the control matrix is unknown, with B = NULL.

Parameters
[in]kThe number of modes to keep after doing SVD.
[in]B(Optional) The control matrix B.

Definition at line 225 of file DMDc.cpp.

Friends And Related Function Documentation

◆ getParametricDMDc

void getParametricDMDc ( DMDc *&  parametric_dmdc,
std::vector< Vector * > &  parameter_points,
std::vector< DMDc * > &  dmdcs,
std::vector< Matrix * >  controls,
Matrix *&  controls_interpolated,
Vector desired_point,
std::string  rbf,
std::string  interp_method,
double  closest_rbf_val,
bool  reorthogonalize_W 
)
friend

Obtain DMD model interpolant at desired parameter point by interpolation of DMD models from training parameter points.

Parameters
[in]parametric_dmdcThe interpolant DMD model at the desired point.
[in]parameter_pointsThe training parameter points.
[in]dmdcsThe DMD objects associated with each training parameter point.
[in]controlsThe matrices of controls from previous runs which we use to interpolate.
[in]controls_interpolatedThe interpolated controls.
[in]desired_pointThe desired point at which to create a parametric DMD.
[in]rbfThe RBF type ("G" == gaussian, "IQ" == inverse quadratic, "IMQ" == inverse multiquadric)
[in]interp_methodThe interpolation method type ("LS" == linear solve, "IDW" == inverse distance weighting, "LP" == lagrangian polynomials)
[in]closest_rbf_valThe RBF parameter determines the width of influence. Set the RBF value of the nearest two parameter points to a value between 0.0 to 1.0
[in]reorthogonalize_WWhether to reorthogonalize the interpolated W (basis) matrix.

The documentation for this class was generated from the following files: