preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
QRFactorization.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <fstream>
5#include <limits>
6#include <string>
7#include <vector>
8#include "logging/Logger.hpp"
10
12
22public:
28 int filter = 0,
29 double omega = 0,
30 double theta = 1. / 0.7,
31 double sigma = std::numeric_limits<double>::min());
32
38 Eigen::MatrixXd A,
39 int filter,
40 double omega = 0,
41 double theta = 1. / 0.7,
42 double sigma = std::numeric_limits<double>::min());
43
49 Eigen::MatrixXd Q,
50 Eigen::MatrixXd R,
51 int rows,
52 int cols,
53 int filter,
54 double omega = 0,
55 double theta = 1. / 0.7,
56 double sigma = std::numeric_limits<double>::min());
57
61 virtual ~QRFactorization() = default;
62
68 void resetFilter(double singularityLimit, std::vector<int> &delIndices, Eigen::MatrixXd &V);
69
73 void reset();
74
78 void reset(
79 Eigen::MatrixXd const &Q,
80 Eigen::MatrixXd const &R,
81 int rows,
82 int cols,
83 double omega = 0,
84 double theta = 1. / 0.7,
85 double sigma = std::numeric_limits<double>::min());
86
90 void reset(
91 Eigen::MatrixXd const &A,
92 int globalRows,
93 double omega = 0,
94 double theta = 1. / 0.7,
95 double sigma = std::numeric_limits<double>::min());
96
101 bool insertColumn(int k, const Eigen::VectorXd &v, double singularityLimit = 0);
102
107 void deleteColumn(int k);
108
114 void pushFront(const Eigen::VectorXd &v);
115
121 void pushBack(const Eigen::VectorXd &v);
122
127 void popFront();
128
133 void popBack();
134
140 void applyFilter(double singularityLimit, std::vector<int> &delIndices, Eigen::MatrixXd &V);
141
145 Eigen::MatrixXd &matrixQ();
146
150 Eigen::MatrixXd &matrixR();
151
152 // @brief returns the number of columns in the QR-decomposition
153 int cols() const;
154 // @brief returns the number of rows in the QR-decomposition
155 int rows() const;
156
157 // @brief optional file-stream for logging output
158 void setfstream(std::fstream *stream);
159
160 // @brief set number of global rows
161 void setGlobalRows(int gr);
162
163 // @brief sets the filtering technique to maintain good conditioning of the least squares system
164 void setFilter(int filter);
165
166 // @brief returns the number of times the QR2 filter step was performed
167 size_t getResetFilterCounter() const;
168
169private:
170 struct givensRot {
171 int i, j;
172 double sigma, gamma;
173 };
174
186 int orthogonalize_stable(Eigen::VectorXd &v, Eigen::VectorXd &r, double &rho, int colNum);
187
200 int orthogonalize(Eigen::VectorXd &v, Eigen::VectorXd &r, double &rho, int colNum);
201
205 void computeReflector(givensRot &grot, double &x, double &y);
206
211 void applyReflector(const givensRot &grot, int k, int l, Eigen::VectorXd &p, Eigen::VectorXd &q);
212
213 logging::Logger _log{"acceleration::QRFactorization"};
214
215 Eigen::MatrixXd _Q;
216 Eigen::MatrixXd _R;
217
218 int _rows;
219 int _cols;
220
222 double _omega;
223 double _theta;
224 double _sigma;
225
226 bool _computeQR2Filter = true; // flag to indicate if the QR2 filter step should be performed
227 size_t _resetFilterCounter = 0; // counter for the number of times the QR2 filter step was performed
228
229 // @brief optional infostream that writes information to file
232
234};
235
236} // namespace precice::acceleration::impl
Class that provides functionality for a dynamic QR-decomposition, that can be updated in O(mn) flops ...
int orthogonalize_stable(Eigen::VectorXd &v, Eigen::VectorXd &r, double &rho, int colNum)
assuming Q(1:n,1:m) has nearly orthonormal columns, this procedure orthogonlizes v(1:n) to the column...
bool insertColumn(int k, const Eigen::VectorXd &v, double singularityLimit=0)
inserts a new column at arbitrary position and updates the QR factorization This function works on th...
int orthogonalize(Eigen::VectorXd &v, Eigen::VectorXd &r, double &rho, int colNum)
assuming Q(1:n,1:m) has nearly orthonormal columns, this procedure orthogonlizes v(1:n) to the column...
Eigen::MatrixXd & matrixR()
returns a matrix representation of the upper triangular matrix R
void applyFilter(double singularityLimit, std::vector< int > &delIndices, Eigen::MatrixXd &V)
filters the least squares system, i.e., the decomposition Q*R = V according to the defined filter tec...
void computeReflector(givensRot &grot, double &x, double &y)
computes parameters for givens matrix G for which (x,y)G = (z,0). replaces (x,y) by (z,...
void resetFilter(double singularityLimit, std::vector< int > &delIndices, Eigen::MatrixXd &V)
Performs a reset of A = QR using the QR2 Filter. This eliminates columns during the reconstruction of...
void pushBack(const Eigen::VectorXd &v)
inserts a new column at position _cols-1, i.e., appends a column at the end and updates the QR factor...
void reset()
resets the QR factorization zo zero Q(0:0, 0:0)R(0:0, 0:0)
void popBack()
deletes the column at position _cols-1, i.e., deletes the last column and updates the QR factorizatio...
void pushFront(const Eigen::VectorXd &v)
inserts a new column at position 0, i.e., shifts right and inserts at first position and updates the ...
virtual ~QRFactorization()=default
Destructor, empty.
QRFactorization(int filter=0, double omega=0, double theta=1./0.7, double sigma=std::numeric_limits< double >::min())
Constructor.
void applyReflector(const givensRot &grot, int k, int l, Eigen::VectorXd &p, Eigen::VectorXd &q)
this procedure replaces the two column matrix [p(k:l-1), q(k:l-1)] by [p(k:l), q(k:l)]*G,...
Eigen::MatrixXd & matrixQ()
returns a matrix representation of the orthogonal matrix Q
void popFront()
deletes the column at position 0, i.e., deletes and shifts columns to the left and updates the QR fac...
void deleteColumn(int k)
updates the factorization A=Q[1:n,1:m]R[1:m,1:n] when the kth column of A is deleted....
This class provides a lightweight logger.
Definition Logger.hpp:17