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 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150
|
// Example reported in Issue #115
// https://github.com/yixuan/spectra/issues/115
#include <iostream>
#include <Eigen/Core>
#include <Eigen/SparseCore>
#include <Spectra/SymGEigsSolver.h>
#include <Spectra/MatOp/SparseCholesky.h>
#include "catch.hpp"
using namespace Spectra;
using Matrix = Eigen::MatrixXd;
using Vector = Eigen::VectorXd;
using SpMat = Eigen::SparseMatrix<double>;
using TriVec = std::vector<Eigen::Triplet<double>>;
using OpType = SparseSymMatProd<double>;
using BOpType = SparseCholesky<double>;
// (Highest) Eigenvalues of: ([A] + lam [B])x = 0
// A := M, (pos. semi definite)
// B := [C + shift*M] (pos. definite)
// -> f = sqrt(1/lam - s)/2pi
void run_test(SpMat& M, SpMat& C, double shift, int nef, int ncv)
{
SpMat A = M;
SpMat B = C + shift * M;
INFO("A =\n " << Matrix(A));
INFO("B =\n " << Matrix(B));
OpType op(A);
BOpType Bop(B);
REQUIRE(Bop.info() == CompInfo::Successful);
SymGEigsSolver<OpType, BOpType, GEigsMode::Cholesky> eigs(op, Bop, nef, ncv);
eigs.init();
int nconv = eigs.compute(SortRule::LargestMagn);
int niter = eigs.num_iterations();
int nops = eigs.num_operations();
INFO("nconv = " << nconv);
INFO("niter = " << niter);
INFO("nops = " << nops);
REQUIRE(eigs.info() == CompInfo::Successful);
Vector evals = eigs.eigenvalues();
Vector lam = 1 / (eigs.eigenvalues().array()) - shift;
Matrix evecs = eigs.eigenvectors();
Matrix resid = A.template selfadjointView<Eigen::Lower>() * evecs -
B.template selfadjointView<Eigen::Lower>() * evecs * evals.asDiagonal();
const double err = resid.array().abs().maxCoeff();
INFO("evals = " << evals.transpose());
INFO("lam = " << lam.transpose());
INFO("U'BU =\n " << evecs.transpose() * B * evecs);
INFO("||AU - BUD||_inf = " << err);
REQUIRE(err == Approx(0.0).margin(1e-9));
}
TEST_CASE("Example #115, case 1", "[example_115]")
{
TriVec C_tri = {
{0, 0, 1.1807575e+08},
{1, 1, 304744.5},
{1, 5, -152372.25},
{2, 2, 304744.5},
{2, 4, 152372.25},
{3, 3, 15403.85},
{4, 2, 152372.25},
{4, 4, 101581.5},
{5, 1, -152372.25},
{5, 5, 101581.5},
};
TriVec M_tri = {
{0, 0, 1000.0},
{1, 1, 1000.0},
{2, 2, 1000.0},
};
SpMat C(6, 6);
C.setFromTriplets(C_tri.begin(), C_tri.end());
SpMat M(6, 6);
M.setFromTriplets(M_tri.begin(), M_tri.end());
const double shift = 1.0e5;
run_test(M, C, shift, 4, 5);
}
TEST_CASE("Example #115, case 2", "[example_115]")
{
TriVec C_tri = {
{0, 0, 2.361515e+08},
{1, 1, 609489.01},
{2, 2, 609489.01},
{3, 3, 30807.7},
{4, 4, 203163},
{5, 5, 203163},
{6, 0, -1.1807575e+08},
{0, 6, -1.1807575e+08},
{6, 6, 1.1807575e+08},
{7, 1, -304744.5},
{1, 7, -304744.5},
{7, 5, -152372.25},
{5, 7, -152372.25},
{7, 7, 304744.5},
{8, 2, -304744.5},
{2, 8, -304744.5},
{8, 4, 152372.25},
{4, 8, 152372.25},
{8, 8, 304744.5},
{9, 3, -15403.85},
{3, 9, -15403.85},
{9, 9, 15403.85},
{10, 2, -152372.25},
{2, 10, -152372.25},
{10, 4, 50790.751},
{4, 10, 50790.751},
{10, 8, 152372.25},
{8, 10, 152372.25},
{10, 10, 101581.5},
{11, 1, 152372.25},
{1, 11, 152372.25},
{11, 5, 50790.751},
{5, 11, 50790.751},
{11, 7, -152372.25},
{7, 11, -152372.25},
{11, 11, 101581.5},
};
TriVec M_tri = {
{0, 0, 1000.0},
{1, 1, 1000.0},
{2, 2, 1000.0},
};
SpMat C(12, 12);
C.setFromTriplets(C_tri.begin(), C_tri.end());
SpMat M(12, 12);
M.setFromTriplets(M_tri.begin(), M_tri.end());
const double shift = 1.0e5;
run_test(M, C, shift, 5, 8);
}
|