OMDc
1.0.0
created from 0b26fa0 on deploy-documentation
|
interface to internal representation of OMDc data More...
#include <DataStore.hpp>
Public Types | |
typedef Eigen::Block< const Eigen::MatrixXd, Eigen::Dynamic, Eigen::Dynamic, true > | Block |
short-hand for matrix block expression | |
typedef Eigen::Block< const Eigen::Ref< const Eigen::MatrixXd >, Eigen::Dynamic, Eigen::Dynamic, true > | RefBlock |
short-hand for referenced matrix block | |
Public Member Functions | |
const RefBlock | X () const |
read-access to 'left' snaps More... | |
const RefBlock | Y () const |
read-access to 'right' snaps More... | |
const Block | LTX () const |
read-access to projected input 'left' snaps More... | |
const Block | LTY () const |
read-access to projected output 'right' snaps More... | |
DataStore (const Eigen::Ref< const Eigen::MatrixXd > S, const Eigen::Ref< const Eigen::MatrixXd > U, Eigen::Index r) | |
standard contructor from snapshots and inputs More... | |
void | updateSystemMatrices (const Eigen::Ref< const Eigen::MatrixXd > L) |
recomputes internal representations for given modes More... | |
const Eigen::MatrixXd | getSystemMatrix () const |
system matrix of linear system from internal representation More... | |
const Eigen::MatrixXd & | getProjector () const |
get state projector More... | |
const Eigen::MatrixXd | getInputMatrix () const |
input matrix of linear system from internal representation More... | |
const double | getResidual () const |
residual of squared Fobenius norm || Y - LML'X - LPU || More... | |
const Eigen::VectorXd | getInitialState () const |
computes intial state from identified system More... | |
Public Attributes | |
const Eigen::Index | m |
number of snapshots | |
const Eigen::Index | r |
rank of reduced space | |
const Eigen::Index | n |
number of spatial locations (n > m) | |
const Eigen::Index | p |
number of inputs | |
interface to internal representation of OMDc data
The original snapshot data is held by reference. This class actually stores only three projection matrices:
The projectors are obtained from QR decompositions, and thus, ensure numerical accuracy by avoiding summation errors.
OMDc::DataStore::DataStore | ( | const Eigen::Ref< const Eigen::MatrixXd > | S, |
const Eigen::Ref< const Eigen::MatrixXd > | U, | ||
Eigen::Index | r | ||
) |
standard contructor from snapshots and inputs
S | snapshots [n x m] |
U | inputs [p x (m-1)] |
r | dimension of reduced system |
const Eigen::VectorXd OMDc::DataStore::getInitialState | ( | ) | const |
computes intial state from identified system
The initial reduced state is computed as the solution to
The initial state is found such that the identified model evolves as closely as possible to the data set. Observe that the OMDc problem states the residual differently: it defines a stepwise error.
Eigen::Vector
const Eigen::MatrixXd OMDc::DataStore::getInputMatrix | ( | ) | const |
input matrix of linear system from internal representation
The input matrix is computed by
using SVD for the pseudo-inverse
const Eigen::MatrixXd& OMDc::DataStore::getProjector | ( | ) | const |
get state projector
The projector is represented by which is defined as
const double OMDc::DataStore::getResidual | ( | ) | const |
residual of squared Fobenius norm || Y - LML'X - LPU ||
The cost function of OMDc is simplified to
using the projector stored internally.
const Eigen::MatrixXd OMDc::DataStore::getSystemMatrix | ( | ) | const |
system matrix of linear system from internal representation
The system matrix is computed by
using SVD for the pseudo-inverse
|
inline |
read-access to projected input 'left' snaps
Eigen::Block
L'X
|
inline |
read-access to projected output 'right' snaps
Eigen::Block
L'Y void OMDc::DataStore::updateSystemMatrices | ( | const Eigen::Ref< const Eigen::MatrixXd > | L | ) |
recomputes internal representations for given modes
L | projection matrix [n x r] |
|
inline |
read-access to 'left' snaps
Eigen::Block
to 'left' snaps
|
inline |
read-access to 'right' snaps
Eigen::Block
to 'right' snaps