File: ComplexEigs.cpp

package info (click to toggle)
spectra 1.2.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 1,788 kB
  • sloc: cpp: 23,044; ansic: 175; fortran: 131; makefile: 90
file content (177 lines) | stat: -rw-r--r-- 4,461 bytes parent folder | download
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
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
#include <iostream>
#include <type_traits>
#include <random>
#include <complex>
#include <Eigen/Core>
#include <Eigen/SparseCore>
#include <Spectra/GenEigsSolver.h>
#include <Spectra/MatOp/DenseGenMatProd.h>
#include <Spectra/MatOp/SparseGenMatProd.h>

#include "catch.hpp"

using namespace Spectra;

using Complex = std::complex<double>;
using ComplexMatrix = Eigen::MatrixXcd;
using ComplexVector = Eigen::VectorXcd;
using ComplexSpMatrix = Eigen::SparseMatrix<Complex>;

ComplexSpMatrix gen_sparse_data(int n, double prob = 0.5)
{
    ComplexSpMatrix mat(n, n);
    std::default_random_engine gen;
    gen.seed(0);
    std::uniform_real_distribution<double> distr(0.0, 1.0);
    for (int i = 0; i < n; i++)
    {
        for (int j = 0; j < n; j++)
        {
            if (distr(gen) < prob)
            {
                double re = distr(gen) - 0.5;
                double im = distr(gen) - 0.5;
                mat.insert(i, j) = Complex(re, im);
            }
        }
    }
    return mat;
}

template <typename MatType, typename Solver>
void run_test(const MatType& mat, Solver& eigs, SortRule selection, bool allow_fail = false)
{
    eigs.init();
    // maxit = 300 to reduce running time for failed cases
    int nconv = eigs.compute(selection, 300);
    int niter = eigs.num_iterations();
    int nops = eigs.num_operations();

    if (allow_fail && eigs.info() != CompInfo::Successful)
    {
        WARN("FAILED on this test");
        std::cout << "nconv = " << nconv << std::endl;
        std::cout << "niter = " << niter << std::endl;
        std::cout << "nops  = " << nops << std::endl;
        return;
    }
    else
    {
        INFO("nconv = " << nconv);
        INFO("niter = " << niter);
        INFO("nops  = " << nops);
        REQUIRE(eigs.info() == CompInfo::Successful);
    }

    ComplexVector evals = eigs.eigenvalues();
    ComplexMatrix evecs = eigs.eigenvectors();

    ComplexMatrix resid = mat * evecs - evecs * evals.asDiagonal();
    const double err = resid.array().abs().maxCoeff();

    INFO("||AU - UD||_inf = " << err);
    REQUIRE(err == Approx(0.0).margin(1e-9));
}

template <typename MatType>
void run_test_sets(const MatType& A, int k, int m)
{
    constexpr bool is_dense = std::is_same<MatType, ComplexMatrix>::value;
    using DenseOp = DenseGenMatProd<Complex>;
    using SparseOp = SparseGenMatProd<Complex>;
    using OpType = typename std::conditional<is_dense, DenseOp, SparseOp>::type;

    OpType op(A);
    GenEigsSolver<OpType> eigs(op, k, m);

    SECTION("Largest Magnitude")
    {
        run_test(A, eigs, SortRule::LargestMagn);
    }
    SECTION("Largest Real Part")
    {
        run_test(A, eigs, SortRule::LargestReal);
    }
    SECTION("Largest Imaginary Part")
    {
        run_test(A, eigs, SortRule::LargestImag);
    }
    SECTION("Smallest Magnitude")
    {
        run_test(A, eigs, SortRule::SmallestMagn, true);
    }
    SECTION("Smallest Real Part")
    {
        run_test(A, eigs, SortRule::SmallestReal);
    }
    SECTION("Smallest Imaginary Part")
    {
        run_test(A, eigs, SortRule::SmallestImag, true);
    }
}

TEST_CASE("Eigensolver of general complex matrix [10x10]", "[eigs_complex]")
{
    std::srand(123);

    const ComplexMatrix A = ComplexMatrix::Random(10, 10);
    int k = 3;
    int m = 6;

    run_test_sets(A, k, m);
}

TEST_CASE("Eigensolver of general complex matrix [100x100]", "[eigs_complex]")
{
    std::srand(123);

    const ComplexMatrix A = ComplexMatrix::Random(100, 100);
    int k = 10;
    int m = 30;

    run_test_sets(A, k, m);
}

TEST_CASE("Eigensolver of general complex matrix [1000x1000]", "[eigs_complex]")
{
    std::srand(123);

    const ComplexMatrix A = ComplexMatrix::Random(1000, 1000);
    int k = 20;
    int m = 50;

    run_test_sets(A, k, m);
}

TEST_CASE("Eigensolver of sparse complex matrix [10x10]", "[eigs_complex]")
{
    std::srand(123);

    const ComplexSpMatrix A = gen_sparse_data(10, 0.5);
    int k = 3;
    int m = 6;

    run_test_sets(A, k, m);
}

TEST_CASE("Eigensolver of sparse complex matrix [100x100]", "[eigs_complex]")
{
    std::srand(123);

    const ComplexSpMatrix A = gen_sparse_data(100, 0.1);
    int k = 10;
    int m = 30;

    run_test_sets(A, k, m);
}

TEST_CASE("Eigensolver of sparse complex matrix [1000x1000]", "[eigs_complex]")
{
    std::srand(123);

    const ComplexSpMatrix A = gen_sparse_data(1000, 0.01);
    int k = 20;
    int m = 50;

    run_test_sets(A, k, m);
}