OMDc  2.0.0
created from 2320578 on deploy-documentation
DataStore.hpp
1 // Copyright 2026 Lucas Mieg, Martin Moennigmann
2 //
3 // Licensed under the Apache License, Version 2.0 (the "License");
4 // you may not use this file except in compliance with the License.
5 // You may obtain a copy of the License at
6 //
7 // http://www.apache.org/licenses/LICENSE-2.0
8 //
9 // Unless required by applicable law or agreed to in writing, software
10 // distributed under the License is distributed on an "AS IS" BASIS,
11 // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
12 // See the License for the specific language governing permissions and
13 // limitations under the License.
14 
15 #ifndef __INC_OMD_DATA_STORE__
16 #define __INC_OMD_DATA_STORE__
17 
18 #include <Eigen/Core>
19 
21 namespace OMD {
22 
24 class DataStore
25 {
26 public:
28  const Eigen::Index m;
30  const Eigen::Index r;
32  const Eigen::Index n;
33 protected:
35  const Eigen::Ref<const Eigen::MatrixXd> S;
37  Eigen::MatrixXd LTS;
39  const double normY;
41  Eigen::MatrixXd W;
43  Eigen::MatrixXd XTL_qr_rep;
45  Eigen::MatrixXd YTL_qr_rep_top;
47  Eigen::MatrixXd YTL_qr_rep_bottom;
48 public:
50  typedef Eigen::Block<const Eigen::MatrixXd,Eigen::Dynamic,Eigen::Dynamic,true> Block;
52  typedef Eigen::Block<const Eigen::Ref<const Eigen::MatrixXd>,Eigen::Dynamic,Eigen::Dynamic,true> RefBlock;
55  inline const RefBlock X() const {
56  return S.leftCols(m-1);
57  };
60  inline const RefBlock Y() const {
61  return S.rightCols( m-1);
62  };
65  inline const Block LTX() const {
66  return LTS.leftCols(m-1);
67  };
70  inline const Block LTY() const {
71  return LTS.rightCols(m-1);
72  };
78  DataStore( const Eigen::Ref<const Eigen::MatrixXd> S, Eigen::Index r);
90  void updateSystemMatrices(const Eigen::Ref<const Eigen::MatrixXd> L);
104  const Eigen::MatrixXd getSystemMatrix() const;
116  const double getResidual() const;
132  inline const Eigen::MatrixXd& getProjector() const {return W;};
136  inline const Eigen::MatrixXd& getYstar() const {return YTL_qr_rep_bottom;};
137 };
138 
139 };
140 
141 #endif
OMD::DataStore
interface to internal representation of OMD data
Definition: DataStore.hpp:25
OMD::DataStore::m
const Eigen::Index m
number of snapshots
Definition: DataStore.hpp:28
OMD::DataStore::updateSystemMatrices
void updateSystemMatrices(const Eigen::Ref< const Eigen::MatrixXd > L)
recomputes internal representations for given modes
OMD::DataStore::getYstar
const Eigen::MatrixXd & getYstar() const
part of Y'L that is orthogonal to X'L as R from QR
Definition: DataStore.hpp:136
OMD::DataStore::DataStore
DataStore(const Eigen::Ref< const Eigen::MatrixXd > S, Eigen::Index r)
standard contructor from snapshots
OMD::DataStore::getSystemMatrix
const Eigen::MatrixXd getSystemMatrix() const
system matrix of linear system M = (L'YX'L) inv(L'XX'L)
OMD::DataStore::X
const RefBlock X() const
read-access to left snaps
Definition: DataStore.hpp:55
OMD::DataStore::n
const Eigen::Index n
number of spatial locations (n > m)
Definition: DataStore.hpp:32
OMD::DataStore::LTS
Eigen::MatrixXd LTS
product of current mode and snapshots [r x m]
Definition: DataStore.hpp:37
OMD
Optimal Mode Decomposition.
Definition: DataStore.hpp:21
OMD::DataStore::XTL_qr_rep
Eigen::MatrixXd XTL_qr_rep
QR representation of X'L as [r x r] (upper triangular)
Definition: DataStore.hpp:43
OMD::DataStore::normY
const double normY
Frobenius norm of left matrices.
Definition: DataStore.hpp:39
OMD::DataStore::getResidual
const double getResidual() const
squared Frobenius residual: ||Y|| - ||(L'YX'L) inv(L'XX'L) (L'XY'L)||
OMD::DataStore::RefBlock
Eigen::Block< const Eigen::Ref< const Eigen::MatrixXd >, Eigen::Dynamic, Eigen::Dynamic, true > RefBlock
short-hand for referenced matrix block
Definition: DataStore.hpp:52
OMD::DataStore::Y
const RefBlock Y() const
read-access to right snaps
Definition: DataStore.hpp:60
OMD::DataStore::LTX
const Block LTX() const
read-access to projected left snaps
Definition: DataStore.hpp:65
OMD::DataStore::W
Eigen::MatrixXd W
projector such that I - Z pinv(Z) = W W' with Z = X'L
Definition: DataStore.hpp:41
OMD::DataStore::Block
Eigen::Block< const Eigen::MatrixXd, Eigen::Dynamic, Eigen::Dynamic, true > Block
short-hand for matrix block expression
Definition: DataStore.hpp:50
OMD::DataStore::YTL_qr_rep_bottom
Eigen::MatrixXd YTL_qr_rep_bottom
QR representation of Y'L as [r x r] (upper triangular)
Definition: DataStore.hpp:47
OMD::DataStore::LTY
const Block LTY() const
read-access to projected right snaps
Definition: DataStore.hpp:70
OMD::DataStore::getProjector
const Eigen::MatrixXd & getProjector() const
projector onto orthog. space of X'L
Definition: DataStore.hpp:132
OMD::DataStore::r
const Eigen::Index r
rank of reduced space
Definition: DataStore.hpp:30
OMD::DataStore::S
const Eigen::Ref< const Eigen::MatrixXd > S
reference to snapshots (no data is stored)
Definition: DataStore.hpp:35
OMD::DataStore::YTL_qr_rep_top
Eigen::MatrixXd YTL_qr_rep_top
QR representation of Y'L as [r x r] (top part)
Definition: DataStore.hpp:45
{}