comp: Compose two functions created using either sourceCppAD or...

Description Usage Arguments Value References Examples

Description

Returns a function with a matrix X as input which computes the value of (f \circ g)(x) = f(g(x)). Note that the order of the composition is such that g is applied to X first. f is then applied to the result of g. The returned function is compatible with both J and H, and if J and H are applied to a function produced by composition, the resulting Jacobian or Hessian matrices are constructed from the Jacobians and Hessians of f and g using a combinatorical form of Faa di Bruno's formula (Hardy 2006). The functions f and g must both be functions of a single argument. Note that a function of multiple arguments can be Curried into a function with a single argument. The R package functional provides the method Curry which is convenient for this purpose.

Usage

1
f %.% g

Arguments

f

Function to be composed with g

g

Function to be composed with f

Value

A function which computes the composition f and g

References

\insertRef

Hardy06RcppEigenAD

\insertRef

ma-tsoy-wo-faa-di-brunoRcppEigenAD

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
library(RcppEigenAD)
# define a function to calculate the eigen vectors of a matrix
f<-sourceCppAD('
ADmat f(const ADmat& X)
{
  Eigen::EigenSolver<ADmat> es(X);
  return es.pseudoEigenvectors();

}
')
# define function to calculate the inverse of a matrix
g<-sourceCppAD('
ADmat g(const ADmat& X)
{
  return X.inverse(); 
}
')
# compose f and g to produce a functions to calculate the eigenvectors of the inverse of a matrix
h<-f%.%g  #  h = f o g
x<-matrix(c(1,2,3,4),2,2)
x<-x%*%t(x) # positive definite matrix
J(h)(x) # Jacobian of h = f o g
H(h)(x) # Stacked Hessians of h = f o g
# redefine h as a function to directly calculate the eigenvectors of the inverse of a matrix
h<-sourceCppAD('
ADmat h(const ADmat& X)
{
  Eigen::EigenSolver<ADmat> es(X.inverse());
  return es.pseudoEigenvectors();
}
')
# calculate the Jacobian and Hessian of h to compare with previous result
J(h)(x) # Jacobian of h = f o g
H(h)(x) # Stacked Hessians of h = f o g

RcppEigenAD documentation built on May 2, 2019, 5:34 a.m.