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

#include <AdaptiveDMD.h>

Inheritance diagram for CAROM::AdaptiveDMD:
CAROM::DMD

Public Member Functions

 AdaptiveDMD (int dim, double desired_dt=-1.0, std::string rbf="G", std::string interp_method="LS", double closest_rbf_val=0.9, bool alt_output_basis=false, Vector *state_offset=NULL)
 Constructor. More...
 
 ~AdaptiveDMD ()
 Destroy the AdaptiveDMD object.
 
void train (double energy_fraction, const Matrix *W0=NULL, double linearity_tol=0.0)
 
void train (int k, const Matrix *W0=NULL, double linearity_tol=0.0)
 
double getTrueDt () const
 Get the true dt between interpolated snapshots.
 
const MatrixgetInterpolatedSnapshots ()
 Get the interpolated snapshot matrix contained within d_interp_snapshots.
 
- Public Member Functions inherited from CAROM::DMD
 DMD (int dim, double dt, bool alt_output_basis=false, Vector *state_offset=NULL)
 Constructor. Basic DMD with uniform time step size. More...
 
 DMD (std::string base_file_name)
 Constructor. DMD from saved models. More...
 
virtual ~DMD ()
 Destroy the DMD object.
 
virtual void setOffset (Vector *offset_vector, int order)
 Set the offset of a certain order.
 
virtual void takeSample (double *u_in, double t)
 Sample the new state, u_in. Any samples in d_snapshots taken at the same or later time will be erased. More...
 
void projectInitialCondition (const Vector *init, double t_offset=-1.0)
 Project new initial condition using d_phi. Calculate pinv(phi) x init, or more precisely, (phi* x phi)^{-1} x phi* x init, where phi* is the conjugate transpose. More...
 
Vectorpredict (double t, int deg=0)
 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 DMD record in CSV files.
 

Additional Inherited Members

- Protected Member Functions inherited from CAROM::DMD
 DMD (int dim, bool alt_output_basis=false, Vector *state_offset=NULL)
 Constructor. Variant of DMD with non-uniform time step size. More...
 
 DMD (std::vector< std::complex< double >> eigs, Matrix *phi_real, Matrix *phi_imaginary, int k, double dt, double t_offset, Vector *state_offset)
 Constructor. Specified from DMD components. More...
 
 DMD ()
 Unimplemented default constructor.
 
 DMD (const DMD &other)
 Unimplemented copy constructor.
 
DMDoperator= (const DMD &rhs)
 Unimplemented assignment operator.
 
std::pair< Matrix *, Matrix * > phiMultEigs (double t, int deg=0)
 Internal function to multiply d_phi with the eigenvalues.
 
void constructDMD (const Matrix *f_snapshots, int rank, int num_procs, const Matrix *W0, double linearity_tol)
 Construct the DMD object.
 
virtual std::pair< Matrix *, Matrix * > computeDMDSnapshotPair (const Matrix *snapshots)
 Returns a pair of pointers to the minus and plus snapshot matrices.
 
virtual void computePhi (DMDInternal dmd_internal_obj)
 Compute phi.
 
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, double t=0.0, int deg=0)
 Add the appropriate offset when predicting the solution.
 
const MatrixcreateSnapshotMatrix (std::vector< Vector * > snapshots)
 Get the snapshot matrix contained within d_snapshots.
 
- Protected Attributes inherited from CAROM::DMD
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.
 
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_sampled_times
 The stored times of each sample.
 
Vectord_state_offset = NULL
 State offset in snapshot.
 
bool d_trained
 Whether the DMD 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 DMD modes.
 
int d_k
 The number of columns used after obtaining the SVD decomposition.
 
Matrixd_basis = NULL
 The left singular vector basis.
 
bool d_alt_output_basis = false
 Whether to use the alternative basis for output, i.e. phi = U^(+)*V*Omega^(-1)*X.
 
Matrixd_A_tilde = NULL
 A_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.
 
std::vector< std::complex< double > > d_eigs
 A vector holding the complex eigenvalues of the eigenmodes.
 

Detailed Description

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

Definition at line 31 of file AdaptiveDMD.h.

Constructor & Destructor Documentation

◆ AdaptiveDMD()

CAROM::AdaptiveDMD::AdaptiveDMD ( int  dim,
double  desired_dt = -1.0,
std::string  rbf = "G",
std::string  interp_method = "LS",
double  closest_rbf_val = 0.9,
bool  alt_output_basis = false,
Vector state_offset = NULL 
)

Constructor.

Parameters
[in]dimThe full-order state dimension.
[in]desired_dtThe constant step size for uniform interpolation of samples. If set equal to or below 0.0, desired_dt will be set to the median of the different dt's between the samples.
[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]alt_output_basisWhether to use the alternative basis for output, i.e. phi = U^(+)*V*Omega^(-1)*X.
[in]state_offsetThe state offset.

Definition at line 22 of file AdaptiveDMD.cpp.

Member Function Documentation

◆ train() [1/2]

void CAROM::AdaptiveDMD::train ( double  energy_fraction,
const Matrix W0 = NULL,
double  linearity_tol = 0.0 
)
virtual
Parameters
[in]energy_fractionThe energy fraction to keep after doing SVD.
[in]W0The initial basis to prepend to W.
[in]linearity_tolThe tolerance for determining whether a column of W is linearly independent with W0.

Reimplemented from CAROM::DMD.

Definition at line 47 of file AdaptiveDMD.cpp.

◆ train() [2/2]

void CAROM::AdaptiveDMD::train ( int  k,
const Matrix W0 = NULL,
double  linearity_tol = 0.0 
)
virtual
Parameters
[in]kThe number of modes (eigenvalues) to keep after doing SVD.
[in]W0The initial basis to prepend to W.
[in]linearity_tolThe tolerance for determining whether a column of W is linearly independent with W0.

Reimplemented from CAROM::DMD.

Definition at line 59 of file AdaptiveDMD.cpp.


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