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

interface to internal representation of OMD data More...

#include <DataStore.hpp>

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

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 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 m
 number of snapshots
 
const Eigen::Index r
 rank of reduced space
 
const Eigen::Index n
 number of spatial locations (n > m)
 

Protected Attributes

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 OMD data

Constructor & Destructor Documentation

◆ DataStore()

OMD::DataStore::DataStore ( const Eigen::Ref< const Eigen::MatrixXd >  S,
Eigen::Index  r 
)

standard contructor from snapshots

The internal data member are initialized to correct dimensions.

Parameters
Ssnapshots
rdimension of reduced system

Member Function Documentation

◆ getProjector()

const Eigen::MatrixXd& OMD::DataStore::getProjector ( ) const
inline

projector onto orthog. space of X'L

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

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

This projection is applied to $ Y^\top L $, such that the internal representation renders the result to

\[ W\,W^\top\, Y^\top L = Q \begin{bmatrix} 0 \\ I \end{bmatrix} Y^\ast = \Pi \, Y^\ast \]

where this method gives access to $ \Pi \in \mathbb{R}^{(m-1)\times r} $, see DataStore::updateSystemMatrices() for the definition of the QR decomposition.

Returns
ref. to matrix [(m-1) x r]

◆ getResidual()

const double OMD::DataStore::getResidual ( ) const

squared Frobenius residual: ||Y|| - ||(L'YX'L) inv(L'XX'L) (L'XY'L)||

The cost function of OMD is simplified to

\[ \Vert Y -LML^\top X \Vert_F^2 = \Vert Y \Vert_F^2 - \Vert L^\top Y \Vert_F^2 + \Vert Y^\ast \Vert_F^2 \]

This form may be further simplified to unite the second and third term. Though, the chosen forms makes it reusable for OMDc

Returns
computed residual from current internal state

◆ getSystemMatrix()

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

system matrix of linear system M = (L'YX'L) inv(L'XX'L)

the system matrix computes

\[ M = ( L^\top Y ) \, (L^\top X )^\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()

Returns
system matrix from internal rep. [r x r]

◆ getYstar()

const Eigen::MatrixXd& OMD::DataStore::getYstar ( ) const
inline

part of Y'L that is orthogonal to X'L as R from QR

Note
only essential part of QR decomposition without Q, to be used with DataStore::getProjector()
Returns
ref. to intertal data [r x r]

◆ LTX()

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

read-access to projected left snaps

Returns
Eigen::Block: L'X

◆ LTY()

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

read-access to projected right snaps

Returns
Eigen::Block L'Y

◆ updateSystemMatrices()

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

recomputes internal representations for given modes

The internal representation of the data reads

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

simplifying the computation of OMD::getResidual() and OMD::getSystemMatrix() . The results of the decomposition are stored in protected members.

Parameters
Lprojection matrix

◆ X()

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

read-access to left snaps

Returns
Eigen::Block to left snaps

◆ Y()

const RefBlock OMD::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:
{}