|
OMDc
2.0.0
created from 2320578 on deploy-documentation
|
interface to internal representation of OMDc data More...
#include <DataStore.hpp>
Public Member Functions | |
| 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 | getInputMatrix () const |
| input matrix of linear system from internal representation More... | |
| const Eigen::VectorXd | getInitialState () const |
| computes intial state from identified system More... | |
Public Member Functions inherited from OMD::DataStore | |
| 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 left snaps More... | |
| const Block | LTY () const |
| read-access to projected right snaps More... | |
| DataStore (const Eigen::Ref< const Eigen::MatrixXd > S, Eigen::Index r) | |
| standard contructor from snapshots 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 M = (L'YX'L) inv(L'XX'L) More... | |
| const double | getResidual () const |
| squared Frobenius residual: ||Y|| - ||(L'YX'L) inv(L'XX'L) (L'XY'L)|| More... | |
| const Eigen::MatrixXd & | getProjector () const |
| projector onto orthog. space of X'L More... | |
| const Eigen::MatrixXd & | getYstar () const |
| part of Y'L that is orthogonal to X'L as R from QR More... | |
Public Attributes | |
| const Eigen::Index | p |
| number of inputs | |
Public Attributes inherited from OMD::DataStore | |
| 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) | |
Additional Inherited Members | |
Public Types inherited from OMD::DataStore | |
| 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 | |
Protected Attributes inherited from OMD::DataStore | |
| const Eigen::Ref< const Eigen::MatrixXd > | S |
| reference to snapshots (no data is stored) | |
| Eigen::MatrixXd | LTS |
| product of current mode and snapshots [r x m] | |
| const double | normY |
| Frobenius norm of left matrices. | |
| Eigen::MatrixXd | W |
| projector such that I - Z pinv(Z) = W W' with Z = X'L | |
| Eigen::MatrixXd | XTL_qr_rep |
| QR representation of X'L as [r x r] (upper triangular) | |
| Eigen::MatrixXd | YTL_qr_rep_top |
| QR representation of Y'L as [r x r] (top part) | |
| Eigen::MatrixXd | YTL_qr_rep_bottom |
| QR representation of Y'L as [r x r] (upper triangular) | |
interface to internal representation of OMDc data
This class expands the OMD::DataStore by
It is intended to be used with OMD::LineSearchFunctor . The OMDc cost is simplified to
where the
denotes matrix representation of the projector onto the orthogonal space of the inputs. This representation is identical to the OMD cost but for the last term which is augmented by
. Hence, the internal data is updated by
such that the OMD::DataStore methods can be applied. E.g. the system matrix computes
using SVD for the pseudo-inverse. It is internally simplified to
by the QR decomposition explained in DataStore::updateSystemMatrices()
| 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
| void OMDc::DataStore::updateSystemMatrices | ( | const Eigen::Ref< const Eigen::MatrixXd > | L | ) |
recomputes internal representations for given modes
The internal representation of the data reads
where the dimensions of the the QR decomposition are identical to OMD . Hence, the internal data stores are overwritten to contain the projection onto the space normal to the inputs.
| L | projection matrix [n x r] |
{}