1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99
|
#include <iostream>
#include <Eigen/Core>
#include <Spectra/LinAlg/Orthogonalization.h>
#include "catch.hpp"
using namespace Spectra;
using Eigen::MatrixXd;
using Eigen::VectorXd;
using Eigen::Index;
template <typename Matrix>
void check_orthogonality(const Matrix& basis)
{
const double tol = 1e-12;
Matrix xs = basis.transpose() * basis;
INFO("The orthonormalized basis must fulfill that basis.T * basis = I");
INFO("Matrix is\n " << basis);
INFO("Overlap is\n " << xs);
CHECK(xs.isIdentity(tol));
}
TEST_CASE("complete orthonormalization", "[orthogonalisation]")
{
std::srand(123);
const Index n = 20;
MatrixXd mat = MatrixXd::Random(n, n);
SECTION("MGS")
{
MGS_orthogonalisation(mat);
check_orthogonality(mat);
}
SECTION("GS")
{
GS_orthogonalisation(mat);
check_orthogonality(mat);
}
SECTION("QR")
{
QR_orthogonalisation(mat);
check_orthogonality(mat);
}
SECTION("twice_is_enough")
{
twice_is_enough_orthogonalisation(mat);
check_orthogonality(mat);
}
SECTION("JensWehner")
{
JensWehner_orthogonalisation(mat);
check_orthogonality(mat);
}
}
TEST_CASE("Partial orthonormalization", "[orthogonalisation]")
{
std::srand(123);
const Index n = 20;
const Index sub = 5;
Index start = n - sub;
// Create a n x 20 orthonormal basis
MatrixXd mat = MatrixXd::Random(n, start);
QR_orthogonalisation(mat);
mat.conservativeResize(Eigen::NoChange, n);
mat.rightCols(sub) = MatrixXd::Random(n, sub);
SECTION("MGS")
{
MGS_orthogonalisation(mat, start);
check_orthogonality(mat);
}
SECTION("GS")
{
GS_orthogonalisation(mat, start);
check_orthogonality(mat);
}
SECTION("twice_is_enough")
{
twice_is_enough_orthogonalisation(mat, start);
check_orthogonality(mat);
}
SECTION("JensWehner")
{
JensWehner_orthogonalisation(mat, start);
check_orthogonality(mat);
}
}
|