From 6d48942d7bdb2a015f2a52e7a290a7a9efb17eec Mon Sep 17 00:00:00 2001 From: Jaron Kent-Dobias Date: Wed, 6 Jan 2021 18:17:58 +0100 Subject: Cleaned up code and removed unused portions. --- stereographic.hpp | 95 ++++++++++++++----------------------------------------- 1 file changed, 23 insertions(+), 72 deletions(-) (limited to 'stereographic.hpp') diff --git a/stereographic.hpp b/stereographic.hpp index a3f2cc6..8313f25 100644 --- a/stereographic.hpp +++ b/stereographic.hpp @@ -1,22 +1,21 @@ #include -#include "Eigen/src/Core/util/Meta.h" #include "p-spin.hpp" -#include "unsupported/Eigen/CXX11/src/Tensor/TensorMeta.h" Vector stereographicToEuclidean(const Vector& ζ) { unsigned N = ζ.size() + 1; Vector z(N); + Scalar a = ζ.transpose() * ζ; Scalar b = 2 * sqrt(N) / (1.0 + a); for (unsigned i = 0; i < N - 1; i++) { - z(i) = b * ζ(i); + z(i) = ζ(i); } - z(N - 1) = b * (a - 1.0) / 2.0; + z(N - 1) = (a - 1.0) / 2.0; - return z; + return b * z; } Vector euclideanToStereographic(const Vector& z) { @@ -24,94 +23,46 @@ Vector euclideanToStereographic(const Vector& z) { Vector ζ(N - 1); for (unsigned i = 0; i < N - 1; i++) { - ζ(i) = z(i) / (sqrt(N) - z(N - 1)); + ζ(i) = z(i); } - return ζ; + return ζ / (sqrt(N) - z(N - 1)); } Matrix stereographicJacobian(const Vector& ζ) { - unsigned N = ζ.size() + 1; - Matrix J(N - 1, N); - - Scalar b = ζ.transpose() * ζ; - - for (unsigned i = 0; i < N - 1; i++) { - for (unsigned j = 0; j < N - 1; j++) { - J(i, j) = - 4 * sqrt(N) * ζ(i) * ζ(j) / pow(1.0 + b, 2); - if (i == j) { - J(i, j) += 2 * sqrt(N) * (1.0 + b) / pow(1.0 + b, 2); - } - } - - J(i, N - 1) = 4.0 * sqrt(N) * ζ(i) / pow(1.0 + b, 2); - } - - return J; -} - -Matrix stereographicMetric(const Matrix& J) { - return J * J.adjoint(); -} - -// Gives the Christoffel symbol, with Γ_(i1, i2)^(i3). -Eigen::Tensor stereographicChristoffel(const Vector& ζ, const Matrix& gInvJacConj) { - unsigned N = ζ.size() + 1; - Eigen::Tensor dJ(N - 1, N - 1, N); + unsigned N = ζ.size(); + Matrix J(N, N + 1); Scalar b = 1.0 + (Scalar)(ζ.transpose() * ζ); - for (unsigned i = 0; i < N - 1; i++) { - for (unsigned j = 0; j < N - 1; j++) { - for (unsigned k = 0; k < N - 1; k++) { - dJ(i, j, k) = 16 * sqrt(N) * ζ(i) * ζ(j) * ζ(k) / pow(b, 3); - if (i == j) { - dJ(i, j, k) -= 4 * sqrt(N) * ζ(k) / pow(b, 2); - } - if (i == k) { - dJ(i, j, k) -= 4 * sqrt(N) * ζ(j) / pow(b, 2); - } - if (j == k) { - dJ(i, j, k) -= 4 * sqrt(N) * ζ(i) / pow(b, 2); - } - } - dJ(i, j, N - 1) = - 16 * sqrt(N) * ζ(i) * ζ(j) / pow(b, 3); + for (unsigned i = 0; i < N; i++) { + for (unsigned j = 0; j < N; j++) { + J(i, j) = - ζ(i) * ζ(j); + if (i == j) { - dJ(i, j, N - 1) += 4 * sqrt(N) * ζ(i) / pow(b, 2); + J(i, j) += 0.5 * b; } } - } - std::array, 1> ip = {Eigen::IndexPair(2, 1)}; + J(i, N) = ζ(i); + } - return dJ.contract(Eigen::TensorMap>(gInvJacConj.data(), N - 1, N), ip); + return 4 * sqrt(N + 1) * J / pow(b, 2); } -std::tuple stereographicHamGradHess(const Tensor& J, const Vector& ζ) { - Vector grad; - Matrix hess; - - Matrix jacobian = stereographicJacobian(ζ); - Vector z = stereographicToEuclidean(ζ); - - std::complex hamiltonian; +std::tuple stereographicHamGradHess(const Tensor& J, const Vector& ζ, const Vector& z) { + Scalar hamiltonian; Vector gradZ; Matrix hessZ; std::tie(hamiltonian, gradZ, hessZ) = hamGradHess(J, z); - Matrix metric = stereographicMetric(jacobian); - - grad = metric.llt().solve(jacobian) * gradZ; + Matrix jacobian = stereographicJacobian(ζ); - /* This is much slower to calculate than the marginal speedup it offers... - Eigen::Tensor Γ = stereographicChristoffel(ζ, gInvJac.conjugate()); - Vector df = jacobian * gradZ; - Eigen::Tensor dfT = Eigen::TensorMap>(df.data(), {df.size()}); - std::array, 1> ip = {Eigen::IndexPair(2, 0)}; - Eigen::Tensor H2 = Γ.contract(dfT, ip); - */ + Matrix metric = jacobian * jacobian.adjoint(); - hess = jacobian * hessZ * jacobian.transpose(); // - Eigen::Map(H2.data(), ζ.size(), ζ.size()); + // The metric is Hermitian and positive definite, so a Cholesky decomposition can be used. + Vector grad = metric.llt().solve(jacobian) * gradZ; + Matrix hess = jacobian * hessZ * jacobian.transpose(); return {hamiltonian, grad, hess}; } -- cgit v1.2.3-54-g00ecf