manifold.optim: Manifold optimization

Description Usage Arguments Details Value References Examples

View source: R/manifold_optim.R

Description

Optimize a function on a manifold.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
manifold.optim(
  prob,
  mani.defn,
  method = "LRBFGS",
  mani.params = get.manifold.params(),
  solver.params = get.solver.params(),
  deriv.params = get.deriv.params(),
  x0 = NULL,
  H0 = NULL,
  has.hhr = FALSE
)

moptim(
  prob,
  mani.defn,
  method = "LRBFGS",
  mani.params = get.manifold.params(),
  solver.params = get.solver.params(),
  deriv.params = get.deriv.params(),
  x0 = NULL,
  H0 = NULL,
  has.hhr = FALSE
)

Arguments

prob

Problem definition

mani.defn

Either a Product manifold definition or one of the Manifold definitions

method

Name of optimization method. Currently supported methods are:

  • "LRBFGS": Limited-memory RBFGS

  • "LRTRSR1": Limited-memory RTRSR1

  • "RBFGS": Riemannian BFGS

  • "RBroydenFamily": Riemannian Broyden family

  • "RCG": Riemannian conjugate gradients

  • "RNewton": Riemannian line-search Newton

  • "RSD": Riemannian steepest descent

  • "RTRNewton": Riemannian trust-region Newton

  • "RTRSD": Riemannian trust-region steepest descent

  • "RTRSR1": Riemannian trust-region symmetric rank-one update

  • "RWRBFGS": Riemannian BFGS

See Huang et al (2016a, 2016b) for details.

mani.params

Arguments to configure the manifold. Construct with get.manifold.params

solver.params

Arguments to configure the solver. Construct with get.solver.params

deriv.params

Arguments to configure numerical differentiation for gradient and Hessian, which are used if those functions are not specified. Construct with get.deriv.params

x0

Starting point for optimization. A numeric vector whose dimension matches the total dimension of the overall problem

H0

Initial value of Hessian. A d \times d matrix, where d is the dimension of x0

has.hhr

Indicates whether to apply the idea in Huang et al (2015) section 4.1 (TRUE or FALSE)

Details

moptim is an alias for manifold.optim.

Value

xopt

Point returned by the solver

fval

Value of the function evaluated at xopt

normgf

Norm of the final gradient

normgfgf0

Norm of the gradient at final iterate divided by norm of the gradient at initiate iterate

iter

Number of iterations taken by the solver

num.obj.eval

Number of function evaluations

num.grad.eval

Number of gradient evaluations

nR

Number of retraction evaluations

nV

Number of occasions in which vector transport is first computed

nVp

Number of remaining computations of vector transport (excluding count in nV)

nH

Number of actions of Hessian

elapsed

Elapsed time for the solver (in seconds)

References

Wen Huang, P.A. Absil, K.A. Gallivan, Paul Hand (2016a). "ROPTLIB: an object-oriented C++ library for optimization on Riemannian manifolds." Technical Report FSU16-14, Florida State University.

Wen Huang, Kyle A. Gallivan, and P.A. Absil (2016b). Riemannian Manifold Optimization Library. URL https://www.math.fsu.edu/~whuang2/pdf/USER_MANUAL_for_2016-04-29.pdf

Wen Huang, K.A. Gallivan, and P.A. Absil (2015). A Broyden Class of Quasi-Newton Methods for Riemannian Optimization. SIAM Journal on Optimization, 25(3):1660-1685.

S. Martin, A. Raim, W. Huang, and K. Adragni (2020). "ManifoldOptim: An R Interface to the ROPTLIB Library for Riemannian Manifold Optimization." Journal of Statistical Software, 93(1):1-32.

Examples

 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
## Not run: 
# ----- Example with objective and gradient written in R -----
set.seed(1234)

p <- 5; n <- 150
B <- matrix(rnorm(n*n), nrow=n)
B <- B + t(B)
D <- diag(p:1, p)

tx <- function(x) { matrix(x, n, p) }
f <- function(x) { X <- tx(x); Trace( t(X) %*% B %*% X %*% D ) }
g <- function(x) { X <- tx(x); 2 * B %*% X %*% D }

mod <- Module("ManifoldOptim_module", PACKAGE = "ManifoldOptim")
prob <- new(mod$RProblem, f, g)

x0 <- as.numeric(orthonorm(matrix(rnorm(n*p), nrow=n, ncol=p)))
mani.params <- get.manifold.params(IsCheckParams = TRUE)
solver.params <- get.solver.params(IsCheckParams = TRUE)
mani.defn <- get.stiefel.defn(n, p)

res <- manifold.optim(prob, mani.defn, method = "RTRSR1",
	mani.params = mani.params, solver.params = solver.params, x0 = x0)
print(res)
head(tx(res$xopt))

## End(Not run)
## Not run: 
library(ManifoldOptim)
library(RcppArmadillo)

# ----- Example with objective and gradient written in C++ -----
set.seed(1234)

p <- 5; n <- 150
B <- matrix(rnorm(n*n), nrow=n)
B <- B + t(B) # force symmetric
D <- diag(p:1, p)

# The Problem class is written in C++. Get a handle to it and set it up from R
Rcpp::sourceCpp(code = '
//[[Rcpp::depends(RcppArmadillo,ManifoldOptim)]]
#include <RcppArmadillo.h>
#include <ManifoldOptim.h>

using namespace Rcpp;
using namespace arma;

class BrockettProblem : public MatrixManifoldOptimProblem
{
public:
	BrockettProblem(const arma::mat& B, const arma::mat& D)
	: MatrixManifoldOptimProblem(false, true), m_B(B), m_D(D) { }

	virtual ~BrockettProblem() { }

	double objFun(const arma::mat& X) const {
		return arma::trace(X.t() * m_B * X * m_D);
	}

	arma::mat gradFun(const arma::mat& X) const {
		return 2 * m_B * X * m_D;
	}

	const arma::mat& GetB() const { return m_B; }
	const arma::mat& GetD() const { return m_D; }

private:
	arma::mat m_B;
	arma::mat m_D;
};

RCPP_MODULE(Brockett_module) {
	class_<BrockettProblem>("BrockettProblem")
	.constructor<mat,mat>()
	.method("objFun", &BrockettProblem::objFun)
	.method("gradFun", &BrockettProblem::gradFun)
	.method("GetB", &BrockettProblem::GetB)
	.method("GetD", &BrockettProblem::GetD)
	;
}
')

prob <- new(BrockettProblem, B, D)
X0 <- orthonorm(matrix(rnorm(n*p), nrow=n, ncol=p))
x0 <- as.numeric(X0)
tx <- function(x) { matrix(x, n, p) }
mani.params <- get.manifold.params(IsCheckParams = TRUE)
solver.params <- get.solver.params(DEBUG = 0, Tolerance = 1e-4,
	Max_Iteration = 1000, IsCheckParams = TRUE, IsCheckGradHess = FALSE)
mani.defn <- get.stiefel.defn(n, p)

res <- manifold.optim(prob, mani.defn, method = "RTRSR1",
	mani.params = mani.params, solver.params = solver.params, x0 = x0)
print(res)
head(tx(res$xopt))

## End(Not run)

ManifoldOptim documentation built on Dec. 15, 2021, 1:07 a.m.