OMDc  2.0.0
created from 2320578 on deploy-documentation
OMDc::DataStore Class Reference

interface to internal representation of OMDc data More...

#include <DataStore.hpp>

Inheritance diagram for OMDc::DataStore:
OMD::DataStore

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)
 

Detailed Description

interface to internal representation of OMDc data

This class expands the OMD::DataStore by

  1. the input data held by reference
  2. the orthogonal projector of the inputs
  3. a method for computing the identified input matrix
  4. a method for the initial state

It is intended to be used with OMD::LineSearchFunctor . The OMDc cost is simplified to

\[ \Vert Y - L M L^\top X - LPU \Vert_F^2 = \Vert (I - L L^\top)Y \Vert_F^2 + \Vert L^\top Y Z (I - Z^\top X^\top L (L^\top X Z Z^\top X^\top L)^{-1} L^\top X Z \Vert_F^2 \]

where the $ Z\in\mathbb{R}^{(m-1)\times (m-p-1)} $ 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 $ Z$. Hence, the internal data is updated by $ Z$ such that the OMD::DataStore methods can be applied. E.g. the system matrix computes

\[ M = ( L^\top Y Z ) \, (L^\top X Z )^\dagger \]

using SVD for the pseudo-inverse. It is internally simplified to

\[ M^\top = R_X^\dagger \, R_Y \]

by the QR decomposition explained in DataStore::updateSystemMatrices()

Note
The residual, the system matrix, the state projector, and the the derivatives for conjugate gradient method are identical to OMD .

Constructor & Destructor Documentation

◆ DataStore()

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

Parameters
Ssnapshots [n x m]
Uinputs [p x (m-1)]
rdimension of reduced system

Member Function Documentation

◆ getInitialState()

const Eigen::VectorXd OMDc::DataStore::getInitialState ( ) const

computes intial state from identified system

The initial reduced state $ a_0$ is computed as the solution to

\[ \min_{a_0} \sum_{k=0}^{m-1} \Vert L^\top s_k - a_k \Vert_2^2 \quad \mathrm{s.t.} \,\, a_{k+1} = M a_k + P u_k \]

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.

Note
The optimization problem can be solved by substituting the identified dynamics into the cost function, and then stacking block-wise powers of the system and input matrix. Then, the resulting problem in the initial state may be solved by the pseudo-inverse of the (large) stacked matrices. However, this implementation does not use stacked matrices but operates with small memory foot print directly on the state vectors.
Returns
initial reduced state as Eigen::Vector

◆ getInputMatrix()

const Eigen::MatrixXd OMDc::DataStore::getInputMatrix ( ) const

input matrix of linear system from internal representation

The input matrix is computed by

\[ P = ( L^\top Y - M L^\top X ) \, U^\dagger \]

using SVD for the pseudo-inverse

Returns
input matrix [r x p]

◆ updateSystemMatrices()

void OMDc::DataStore::updateSystemMatrices ( const Eigen::Ref< const Eigen::MatrixXd >  L)

recomputes internal representations for given modes

The internal representation of the data reads

\[ \begin{bmatrix} Z^\top X^\top L & Z^\top Y^\top L \end{bmatrix} = Q \, \begin{bmatrix} R_X & R_Y \\ 0 & Y^\ast \end{bmatrix} \]

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.

Note
The internal data available through OMD::DataStore::getProjector() and OMD::DataStore::getYstar() contains the projected data, while the dimensions are not changed
Parameters
Lprojection matrix [n x r]

The documentation for this class was generated from the following file:
{}