| Title: | Differentiating Through Cone Programs |
|---|---|
| Description: | A port of the 'python' 'diffcp' package. Computes the derivative of the optimal solution map of a convex cone program, treating the program as an implicit function of its data (constraint matrix, offset, objective coefficients, and optionally a quadratic), mirroring Agrawal et al. (2019) <doi:10.48550/arXiv.1904.09043>. |
| Authors: | Balasubramanian Narasimhan [aut, cre], Akshay Agrawal [aut], Shane Barratt [aut], Stephen Boyd [aut], Enzo Busseti [aut], Walaa Moursi [aut] |
| Maintainer: | Balasubramanian Narasimhan <[email protected]> |
| License: | Apache License (>= 2) |
| Version: | 0.1.1 |
| Built: | 2026-06-04 19:31:25 UTC |
| Source: | https://github.com/bnaras/diffcp |
Mirrors diffcp.cones.parse_cone_dict.
parse_cone_dict(cone_dict)parse_cone_dict(cone_dict)
cone_dict |
A named list with keys among |
A list of list(name, size) pairs, in the canonical SCS
cone order (CONES).
Mirrors diffcp.cones.pi.
pi(x, cones, dual = FALSE)pi(x, cones, dual = FALSE)
x |
A numeric vector. |
cones |
A list of |
dual |
If |
A numeric vector of the same length as x.
Mirrors diffcp.cone_program.solve_and_derivative. Solves
minimize c^T x s.t. A x + s = b, s in K
(with optional QP 0.5 x^T P x term) and returns the optimal
(x, y, s) together with closures D (forward) and DT (adjoint)
that map perturbations of (A, b, c, [P]) to perturbations of
(x, y, s) and vice versa.
solve_and_derivative( A, b, c, cone_dict, P = NULL, solve_method = "Clarabel", mode = "lsqr", warm_start = NULL, ... )solve_and_derivative( A, b, c, cone_dict, P = NULL, solve_method = "Clarabel", mode = "lsqr", warm_start = NULL, ... )
A |
A sparse |
b |
A numeric offset vector. |
c |
A numeric objective coefficient vector. |
cone_dict |
A named list with cone sizes (keys among
|
P |
Optional sparse |
solve_method |
One of |
mode |
Differentiation mode: |
warm_start |
Optional warm-start |
... |
Additional control parameters forwarded to the solver. |
A list with elements x, y, s, info, D, DT.
info is the solver-status block returned by the underlying
solver (status, iter, solveTime, pobj, ...); see
solve_only for the same shape.
Mirrors diffcp.cone_program.solve_only. Solves
minimize c^T x (+ 0.5 x^T P x)
subject to A x + s = b, s in K
and returns the optimal (x, y, s). Unlike solve_and_derivative
this function does not build the derivative closures.
solve_only( A, b, c, cone_dict, warm_start = NULL, solve_method = "Clarabel", P = NULL, ... )solve_only( A, b, c, cone_dict, warm_start = NULL, solve_method = "Clarabel", P = NULL, ... )
A |
A sparse |
b |
A numeric offset vector. |
c |
A numeric objective coefficient vector. |
cone_dict |
A named list with cone sizes (keys among
|
warm_start |
Optional warm-start |
solve_method |
One of |
P |
Optional sparse |
... |
Additional control parameters forwarded to the solver. |
A list with elements x, y, s, info.
vec_symm
Inverse of vec_symm
unvec_symm(x, dim)unvec_symm(x, dim)
x |
A numeric vector of length |
dim |
The matrix dimension |
The corresponding n x n symmetric matrix.
Returns the upstream pin metadata for this R port: which
cvxgrp/diffcp version / commit the R sources track, the date of
that commit, and diffcp's own CVXPY runtime constraint. The
authoritative record lives in inst/UPSTREAM.dcf and is read at
call time.
upstream_info()upstream_info()
A named character vector with one element per DCF field
(Upstream, Version, Commit, ShortCommit, Date, URL,
CVXPY, Snapshot, Notes).
upstream_info() upstream_info()[["Version"]]upstream_info() upstream_info()[["Version"]]
Returns a vectorized representation of a symmetric matrix X,
with off-diagonal entries scaled by sqrt(2) to make the SCS-style
dot product <vec_symm(A), vec_symm(B)> = trace(A B) hold.
vec_symm(X)vec_symm(X)
X |
A symmetric matrix. |
Mirrors diffcp.cones.vec_symm (Python).
A numeric vector of length n*(n+1)/2.