ActiveSubspaceMethods
Documentation for ActiveSubspaceMethods.
ActiveSubspaceMethods.ADFunctionWrapper
ActiveSubspaceMethods.ActiveSubspaces
ActiveSubspaceMethods.GaussianizedUniformInputFunction
ActiveSubspaceMethods.ModifiedActiveSubspaces
ActiveSubspaceMethods.ASXX_bound
ActiveSubspaceMethods.ActiveSubspacesXXManopt
ActiveSubspaceMethods.GaussianizedUniformInputFunctionAD
ActiveSubspaceMethods.MCActiveSubspacesInput
ActiveSubspaceMethods.MC_variance
ActiveSubspaceMethods.QuadratureActiveSubspacesInput
ActiveSubspaceMethods.QuadratureActiveSubspacesInput
ActiveSubspaceMethods.create_Ur
ActiveSubspaceMethods.nested_MC_err
ActiveSubspaceMethods.orthogonalize
ActiveSubspaceMethods.ADFunctionWrapper
— MethodActiveSubspaceMethods.ADFunctionWrapper(fcn, d, backend; [init_space])
Wraps the evaluation of a pure Julia function fcn
to use automatic differentiation.
Arguments
fcn(z::AbstractVector)::Float64
function we want to reduce dimension ofd::Int
dimension of inputbackend::ADTypes.AbstractADType
Backend for autodiff (e.g., ReverseDiff, Enzyme, etc.). Must be explicitly installed!init_space::Vector
values to prepare the memory spaces for autodiff, defaultzeros(d)
.
ActiveSubspaceMethods.ActiveSubspaces
— MethodActiveSubspaces(input::AbstractActiveSubspacesInput)
Type to represent traditional Active Subspaces method [Russi, 2010].
ActiveSubspaceMethods.GaussianizedUniformInputFunction
— TypeGaussianizedUniformInputFunction([f0], f1!, bounds; [test_eval])
Take a function that we wish to minimize L2 error w.r.t. uniform measure and transform it to be w.r.t. GaussianizedUniformInputFunction
Arguments
f0(z)::Float64
Returns function eval atz
. Will default to usingf1!
if not providedf1!(grad, z)::Float64
Returns same function eval atz
and puts function gradient atz
intograd
bounds::AbstractVector{<:NTuple{2}}
A sequence of (lower, upper) bounds that is as long as the number of inputs.test_eval::Bool
Whether to test function eval when constructing this functor, defaulttrue
.
ActiveSubspaceMethods.ModifiedActiveSubspaces
— MethodModifiedActiveSubspaces(input::AbstractActiveSubspacesInput)
Type to represent Modified Active Subspaces method [Lee, 2019].
ActiveSubspaceMethods.ASXX_bound
— MethodASXX_bound(as::AbstractActiveSubspacesXX, U_perp::AbstractMatrix)
Evaluate the active subspaces++ bound using the perp subspace.
ActiveSubspaceMethods.ActiveSubspacesXXManopt
— FunctionActiveSubspacesXXManopt(input::AbstractActiveSubspacesInput)
Type to represent Active Subspaces ++, proposed in [Li et al., 2025]. Only usable when Manopt and Manifolds are installed.
ActiveSubspaceMethods.GaussianizedUniformInputFunctionAD
— MethodGaussianizedUniformInputFunctionAD(f0, bounds, backend; [test_eval])
Equivalent to GaussianizedUniformInputFunction(f0, ADFunctionWrapper(f0, d, backend), bounds; test_eval)
ActiveSubspaceMethods.MCActiveSubspacesInput
— MethodMCActiveSubspacesInput(eval_grad_fcn!, d, N; [rand_fcn, rng, corrected])
Construct input for an active subspace method from a function (with inplace gradient) and Monte Carlo integration.
Arguments
eval_grad_fcn!(grad, z)
Returns real-valued function and puts gradient of this function ingrad
d::Int
dimensionN::Int
number of samplesrand_fcn(rng, d, N)
function to sample random numbers usingrng
, defaultrandn
rng::AbstractRNG
Random number generatorcorrected::Bool
Use corrected Monte Carlo estimators downstream, defaulttrue
ActiveSubspaceMethods.MC_variance
— MethodMC_variance(fcn_eval,d,N; [rand_fcn,rng])
Approximate the N
sample variance of function fcn_eval
with d
inputs
ActiveSubspaceMethods.QuadratureActiveSubspacesInput
— MethodQuadratureActiveSubspacesInput(eval_grad_fcn!, d, tensor_order, [quad_fcn1d, verbose])
Construct input for an active subspace method from a function (with inplace gradient) and Monte Carlo integration.
Arguments
eval_grad_fcn!(grad, z)
Returns real-valued function and puts gradient of this function ingrad
d::Int
dimensiontensor_order::Int
number of one-dimensional quadrature points, will givetensor_order^d
quadrature pointsquad_fcn1d(order)::Tuple{Vector,Vector}
One-dimensional quadrature, defaults to normalized Gauss-Hermiteverbose::Bool
Whether to print the number of quadrature points.
ActiveSubspaceMethods.QuadratureActiveSubspacesInput
— MethodQuadratureActiveSubspacesInput(eval_grad_fcn!, d, tensor_order, [quad_fcn1d, verbose])
Construct input for an active subspace method from a function (with inplace gradient) and Monte Carlo integration.
Arguments
eval_grad_fcn!(grad, z)
Returns real-valued function and puts gradient of this function ingrad
tensor_orders::NTuple{d,Int}
number of one-dimensional quadrature points for each dimension, will give number of quad points as product of these ordersquad_fcn1d(order)::Tuple{Vector,Vector}
One-dimensional quadrature, defaults to normalized Gauss-Hermiteverbose::Bool
Whether to print the number of quadrature points.
ActiveSubspaceMethods.create_Ur
— Functioncreate_Ur(U_perp, [initial_guess, rng])
Create a projector $U_r$ that projects onto orthogonal space to given $U_\perp$.
ActiveSubspaceMethods.nested_MC_err
— Methodnested_MC_err(fcn_eval, U_r, U_perp, N_inner, N_outer; [rand_fcn, rng])
Calculates the conditional variance of subspaces using a pick-and-freeze nested MC method.
Arguments
fcn_eval(z::Vector)::Float64
Evaluation of function taking in length-d vectorU_r
a size (d,r) orthogonal matrixU_perp
a size (d,d-r) orthogonal matrixN_inner::Int
number of inner Monte Carlo samplesN_outer::Int
number of outer Monte Carlo samplesrand_fcn(rng,k)::Vector{Float64}
sampler of random numbers, defaultrandn
rng::AbstractRNG
random number generator
ActiveSubspaceMethods.orthogonalize
— Methodorthogonalize(U)
Return matrix with orthonormal columns but same range as U