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

interface to internal representation of OMD 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 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...
 

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)
 

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

Parameters
Ssnapshots
rdimension of reduced system

Member Function Documentation

◆ getProjector()

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

projector onto orthog. space of X'L

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

◆ getResidual()

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

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

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)

Returns
system matrix from internal rep. [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

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:

Automatic Control and Systems Theory, Ruhr-Univeristy Bochum