OMDc  1.0.0
created from 0b26fa0 on deploy-documentation
OMDc::DataStore Class Reference

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
 

Detailed Description

interface to internal representation of OMDc data

The original snapshot data is held by reference. This class actually stores only three projection matrices:

  1. the snapshots projected onto the lifting matrix
  2. the orthogonal projector of the inputs
  3. an additional orthogonal projector for the cost function

The projectors are obtained from QR decompositions, and thus, ensure numerical accuracy by avoiding summation errors.

Note
The matrices of the identified system are only computed on demand by respective methods.

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]

◆ getProjector()

const Eigen::MatrixXd& OMDc::DataStore::getProjector ( ) const

get state projector

The projector is represented by $ZW$ which is defined as

\[ ZW\,W^\top Z^\top = Z \big( I - Z^\top X^\top L (L^\top X Z \, Z^\top X^\top L)^{-1} L^\top X Z \big) Z^\top \qquad Z\,Z^\top = \big(I - U^\top (U U^\top)^{-1} U \big) \]

Returns
read-only reference to internal projection matrix

◆ getResidual()

const double OMDc::DataStore::getResidual ( ) const

residual of squared Fobenius norm || Y - LML'X - LPU ||

The cost function of OMDc is simplified to

\[ \Vert Y -LML^\top X-LPU\Vert_F^2 = \Vert Y \Vert_F^2 - \Vert L^\top Y \Vert_F^2 + \Vert L^\top Y \,Z W \Vert_F^2 \]

using the projector stored internally.

Returns
computed residual from current internal state

◆ getSystemMatrix()

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

system matrix of linear system from internal representation

The system matrix is computed by

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

using SVD for the pseudo-inverse

Returns
system matrix [r x r]

◆ LTX()

const Block OMDc::DataStore::LTX ( ) const
inline

read-access to projected input 'left' snaps

Returns
Eigen::Block L'X

◆ LTY()

const Block OMDc::DataStore::LTY ( ) const
inline

read-access to projected output 'right' snaps

Returns
Eigen::Block L'Y

◆ updateSystemMatrices()

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

recomputes internal representations for given modes

Parameters
Lprojection matrix [n x r]

◆ X()

const RefBlock OMDc::DataStore::X ( ) const
inline

read-access to 'left' snaps

Returns
Eigen::Block to 'left' snaps

◆ Y()

const RefBlock OMDc::DataStore::Y ( ) const
inline

read-access to 'right' snaps

Returns
Eigen::Block to 'right' snaps

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

Automatic Control and Systems Theory, Ruhr-Univeristy Bochum