preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
QRFactorizationTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <cmath>
9#include "testing/Testing.hpp"
10
11BOOST_AUTO_TEST_SUITE(AccelerationTests)
12
13using namespace precice;
14using namespace precice::acceleration;
15using namespace precice::acceleration::impl;
16
18 Eigen::MatrixXd &Q,
19 Eigen::MatrixXd &R,
20 Eigen::MatrixXd &A)
21{
22 Eigen::MatrixXd A_prime = Q * R;
23
24 for (int i = 0; i < A.rows(); i++) {
25 for (int j = 0; j < A.cols(); j++) {
26 BOOST_TEST(testing::equals(A_prime(i, j), A(i, j)));
27 }
28 }
29}
30
31void testQTQequalsIdentity(Eigen::MatrixXd &Q)
32{
33 Eigen::MatrixXd QTQ = Q.transpose() * Q;
34
35 // test if Q^TQ equals identity
36 for (int i = 0; i < QTQ.rows(); i++) {
37 for (int j = 0; j < QTQ.cols(); j++) {
38 if (i == j) {
39 BOOST_TEST(testing::equals(QTQ(i, j), 1.0, 1e-12));
40 } else {
41 BOOST_TEST(testing::equals(QTQ(i, j), 0.0, 1e-12));
42 }
43 }
44 }
45}
46
48BOOST_AUTO_TEST_CASE(testQRFactorization)
49{
51 int m = 6, n = 8;
52 int filter = BaseQNAcceleration::QR1FILTER;
53 Eigen::MatrixXd A(n, m);
54
55 // Set values according to Hilbert matrix.
56 for (int i = 0; i < n; i++) {
57 for (int j = 0; j < m; j++) {
58 double entry = 1.0 / static_cast<double>(i + j + 1);
59 A(i, j) = entry;
60 }
61 }
62
63 // compute QR factorization of A via successive inserting of columns
64 QRFactorization qr_1(A, filter);
65
66 // test if Q^TQ equals identity
68 // test if QR equals A
69 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A);
70
75 // ----------- delete last column ---------------
76 Eigen::MatrixXd A_prime1 = A;
77 Eigen::VectorXd col6 = A.col(m - 1);
78 A_prime1.conservativeResize(n, m - 1);
79 qr_1.deleteColumn(m - 1);
80
81 // test if Q^TQ equals identity
83 // test if QR equals A
84 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A_prime1);
85
86 // ----------- delete first column ---------------
87 Eigen::MatrixXd A_prime2 = A_prime1;
88 Eigen::VectorXd col1 = A.col(0);
89 for (int i = 0; i < A_prime2.rows(); i++)
90 for (int j = 0; j < A_prime2.cols() - 1; j++)
91 A_prime2(i, j) = A_prime2(i, j + 1);
92 A_prime2.conservativeResize(n, A_prime2.cols() - 1);
93 qr_1.deleteColumn(0);
94
95 // test if Q^TQ equals identity
97 // test if QR equals A
98 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A_prime2);
99
100 // ----------- add first column -----------------
101 qr_1.insertColumn(0, col1);
102
103 // test if Q^TQ equals identity
105 // test if QR equals A
106 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A_prime1);
107
108 // ----------- add last column -----------------
109 qr_1.insertColumn(qr_1.cols(), col6); // ??
110 // test if Q^TQ equals identity
112 // test if QR equals A
113 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A);
114
115 // ----------- delete middle column ---------------
116 int k = 3;
117 Eigen::MatrixXd A_prime3 = A;
118 Eigen::VectorXd colk = A.col(k);
119 for (int i = 0; i < A_prime3.rows(); i++)
120 for (int j = k; j < A_prime3.cols() - 1; j++)
121 A_prime3(i, j) = A_prime3(i, j + 1);
122 A_prime3.conservativeResize(n, A_prime3.cols() - 1);
123 qr_1.deleteColumn(k);
124
125 // test if Q^TQ equals identity
127 // test if QR equals A
128 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A_prime3);
129
130 // ----------- add middle column -----------------
131 qr_1.insertColumn(k, colk);
132
133 // test if Q^TQ equals identity
135 // test if QR equals A
136 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A);
137
138 // ------------ reset ----------------------------
139 qr_1.reset();
140 qr_1.reset(A, A.rows());
142 // test if QR equals A
143 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A);
144
145 Eigen::MatrixXd q = qr_1.matrixQ();
146 Eigen::MatrixXd r = qr_1.matrixR();
147 qr_1.reset();
148 qr_1.reset(q, r, q.rows(), q.cols());
150 // test if QR equals A
151 testQRequalsA(qr_1.matrixQ(), qr_1.matrixR(), A);
152}
153
154PRECICE_TEST_SETUP(1_rank)
156{
157 int m = 6, n = 8;
158 int filter2 = BaseQNAcceleration::QR2FILTER;
159 int filter3 = BaseQNAcceleration::QR3FILTER;
160 std::vector<int> delIndices;
161 Eigen::MatrixXd A2(n, m);
162 Eigen::MatrixXd A3(n, m);
163
164 // Set values according to Hilbert matrix.
165 for (int i = 0; i < n; i++) {
166 for (int j = 0; j < m; j++) {
167 double entry = 1.0 / static_cast<double>(i + j + 1);
168 A2(i, j) = entry;
169 A3(i, j) = entry;
170 }
171 }
172
173 // compute QR factorization of A via successive inserting of columns
174 QRFactorization qr_2(A2, filter2);
175 QRFactorization qr_3(A3, filter3);
176
177 BOOST_TEST(precice::testing::equals(qr_2.matrixQ(), qr_3.matrixQ()));
178 qr_2.applyFilter(1e-3, delIndices, A2);
179 qr_3.applyFilter(1e-3, delIndices, A3);
180 BOOST_TEST(precice::testing::equals(qr_2.matrixQ(), qr_3.matrixQ()));
181 BOOST_TEST(qr_2.getResetFilterCounter() == 1);
182 BOOST_TEST(qr_3.getResetFilterCounter() == 1);
183
184 qr_2.applyFilter(1e-3, delIndices, A2);
185 qr_3.applyFilter(1e-3, delIndices, A3);
186 BOOST_TEST(precice::testing::equals(qr_2.matrixQ(), qr_3.matrixQ()));
187 BOOST_TEST(qr_2.getResetFilterCounter() == 2);
188 BOOST_TEST(qr_3.getResetFilterCounter() == 1);
189}
190
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
void testQTQequalsIdentity(Eigen::MatrixXd &Q)
void testQRequalsA(Eigen::MatrixXd &Q, Eigen::MatrixXd &R, Eigen::MatrixXd &A)
BOOST_AUTO_TEST_CASE(testQRFactorization)
#define PRECICE_TEST()
Definition Testing.hpp:39
#define PRECICE_TEST_SETUP(...)
Creates and attaches a TestSetup to a Boost test case.
Definition Testing.hpp:29
Class that provides functionality for a dynamic QR-decomposition, that can be updated in O(mn) flops ...
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...
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 reset()
resets the QR factorization zo zero Q(0:0, 0:0)R(0:0, 0:0)
Eigen::MatrixXd & matrixQ()
returns a matrix representation of the orthogonal matrix Q
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....
contains implementations of acceleration schemes.
boost::test_tools::predicate_result equals(const std::vector< float > &VectorA, const std::vector< float > &VectorB, float tolerance)
equals to be used in tests. Compares two std::vectors using a given tolerance. Prints both operands o...
Definition Testing.cpp:93
Main namespace of the precice library.