summaryrefslogtreecommitdiff
path: root/least_squares.cpp
blob: a4ee7273e49d970ef3d8fbc7954c97bf507d3f23 (plain)
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
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
#include <eigen3/Eigen/Dense>
#include <eigen3/unsupported/Eigen/CXX11/Tensor>
#include <eigen3/unsupported/Eigen/CXX11/TensorSymmetry>
#include <getopt.h>
#include <iomanip>

#include "pcg-cpp/include/pcg_random.hpp"
#include "randutils/randutils.hpp"

using Rng = randutils::random_generator<pcg32>;

using Real = double;
using Vector = Eigen::Matrix<Real, Eigen::Dynamic, 1>;
using Matrix = Eigen::Matrix<Real, Eigen::Dynamic, Eigen::Dynamic>;

class Tensor : public Eigen::Tensor<Real, 3> {
  using Eigen::Tensor<Real, 3>::Tensor;

public:
  Matrix operator*(const Vector& x) const {
    std::array<Eigen::IndexPair<int>, 1> ip20 = {Eigen::IndexPair<int>(2, 0)};
    const Eigen::Tensor<Real, 1> xT = Eigen::TensorMap<const Eigen::Tensor<Real, 1>>(x.data(), x.size());
    const Eigen::Tensor<Real, 2> JxT = contract(xT, ip20);
    return Eigen::Map<const Matrix>(JxT.data(), dimension(0), dimension(1));
  }
};

Matrix operator*(const Eigen::Matrix<Real, 1, Eigen::Dynamic>& x, const Tensor& J) {
  std::array<Eigen::IndexPair<int>, 1> ip00 = {Eigen::IndexPair<int>(0, 0)};
  const Eigen::Tensor<Real, 1> xT = Eigen::TensorMap<const Eigen::Tensor<Real, 1>>(x.data(), x.size());
  const Eigen::Tensor<Real, 2> JxT = J.contract(xT, ip00);
  return Eigen::Map<const Matrix>(JxT.data(), J.dimension(1), J.dimension(2));
}

Vector normalize(const Vector& x) {
  return x * sqrt(x.size() / x.squaredNorm());
}

Vector ∂HFromV∂V(const Vector& V, const Matrix& ∂V) {
  return V.transpose() * ∂V;
}

Vector VFromABJx(const Vector& b, const Matrix& A, const Matrix& Jx, const Vector& x) {
  return b + (A + 0.5 * Jx) * x;
}

class QuadraticModel {
private:
  Tensor J;
  Matrix A;
  Vector b;

  std::tuple<Vector, Matrix, const Tensor&> V_∂V_∂∂V(const Vector& x) const {
    Matrix Jx = J * x;
    Vector V = VFromABJx(b, A, Jx, x);
    Matrix ∂V = A + Jx;
    return {V, ∂V, J};
  }

  std::tuple<Vector, Matrix> ∂H_∂∂H(const Vector& x) const {
    auto [V, ∂V, ∂∂V] = V_∂V_∂∂V(x);
    Vector ∂H = ∂HFromV∂V(V, ∂V);
    Matrix ∂∂H = V.transpose() * ∂∂V + ∂V.transpose() * ∂V;
    return {∂H, ∂∂H};
  }

public:
  unsigned N;
  unsigned M;

  QuadraticModel(unsigned N, unsigned M, Real σ², Real μA, Real μJ, Rng& r) : J(M, N, N), A(M, N), b(M), N(N), M(M) {
    Eigen::StaticSGroup<Eigen::Symmetry<1,2>> sym23;

    for (unsigned k = 0; k < N; k++) {
      for (unsigned j = k; j < N; j++) {
        for (unsigned i = 0; i < M; i++) {
          sym23(J, i, j, k) = r.variate<Real, std::normal_distribution>(0, sqrt(2 * μJ) / N);
        }
      }
    }

    for (Real& Aij : A.reshaped()) {
      Aij = r.variate<Real, std::normal_distribution>(0, sqrt(μA / N));
    }

    for (Real& bi : b) {
      bi = r.variate<Real, std::normal_distribution>(0, sqrt(σ²));
    }
  }

  Real H(const Vector& x) const {
    Vector V = VFromABJx(b, A, J * x, x);
    return 0.5 * V.squaredNorm();
  }

  Vector ∇H(const Vector& x) const {
    auto [V, ∂V, ∂∂V] = V_∂V_∂∂V(x);
    Vector ∂H = ∂HFromV∂V(V, ∂V);
    return ∂H - (∂H.dot(x) / x.squaredNorm()) * x;
  }

  Matrix HessH(const Vector& x) const {
    auto [∂H, ∂∂H] = ∂H_∂∂H(x);
    Matrix P = Matrix::Identity(N, N) - x * x.transpose() / x.squaredNorm();
    return P * ∂∂H * P.transpose() - (x.dot(∂H) / N) * Matrix::Identity(N, N);
  }

  Vector spectrum(const Vector& x) const {
    Matrix hessH  = HessH(x);
    Eigen::EigenSolver<Matrix> eigenS(hessH);
    return eigenS.eigenvalues().real();
  }

  Real λmax(const Vector& x) const {
    return spectrum(x).maxCoeff();
  }
};

Vector gradientAscent(const QuadraticModel& M, const Vector& x₀, Real ε = 1e-13) {
  Vector xₜ = x₀;
  Real α = 1;
  Real Hₜ = M.H(x₀);
  Real m;
  Vector ∇H;

  while (
    ∇H = M.H(xₜ), m = ∇H.squaredNorm(),
    m / M.N > ε
  ) {
    Real Hₜ₊₁;
    Vector xₜ₊₁;

    while(
      xₜ₊₁ = normalize(xₜ + α * ∇H), Hₜ₊₁ = M.H(xₜ₊₁),
      Hₜ₊₁ < Hₜ + 0.5 * α * m
    ) {
      α /= 2;
    }

    xₜ = xₜ₊₁;
    Hₜ = Hₜ₊₁;
    α *= 1.25;
  }

  return xₜ;
}

Vector messagePassing(const QuadraticModel& M, const Vector& x₀) {
  Vector σ = x₀ / x₀.norm();

  for (unsigned i = 0; i < M.N; i++) {
    Vector ∇H = M.H(σ);
    Vector v = ∇H / ∇H.norm();

    σ += v;
  }

  return normalize(σ);
}

int main(int argc, char* argv[]) {
  unsigned N = 10;
  Real α = 1;
  Real σ² = 1;
  Real μA = 1;
  Real μJ = 1;
  unsigned samples = 10;

  int opt;

  while ((opt = getopt(argc, argv, "N:a:s:A:J:n:")) != -1) {
    switch (opt) {
    case 'N':
      N = (unsigned)atof(optarg);
      break;
    case 'a':
      α = atof(optarg);
      break;
    case 's':
      σ² = atof(optarg);
      break;
    case 'A':
      μA = atof(optarg);
      break;
    case 'J':
      μJ = atof(optarg);
      break;
    case 'n':
      samples = atoi(optarg);
      break;
    default:
      exit(1);
    }
  }

  unsigned M = std::round(α * N);

  Rng r;

  Vector x = Vector::Zero(N);
  x(0) = sqrt(N);

  std::cout << std::setprecision(15);

  for (unsigned sample = 0; sample < samples; sample++) {
    QuadraticModel* ls = new QuadraticModel(N, M, σ², μA, μJ, r);
    Vector xGD = gradientAscent(*ls, x);
    std::cout << ls->H(xGD) / N << " " << ls->λmax(xGD) << " ";
    delete ls;

    ls = new QuadraticModel(N, M, σ², μA, μJ, r);
    Vector xMP = messagePassing(*ls, x);
    xMP = gradientAscent(*ls, xMP);
    std::cout << ls->H(xMP) / N << " " << ls->λmax(xMP) << std::endl;
    delete ls;
  }

  return 0;
}