Supported operations
Convex.jl supports the following functions. These functions may be composed according to the DCP composition rules to form new convex, concave, or affine expressions.
*
Base.:*
— MethodBase.:*(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
The binary multiplication operator $x \times y$.
Examples
ulia> x = Variable();
julia> 2 * x
* (affine; real)
├─ [2;;]
└─ real variable (id: 709…007)
julia> x = Variable(3);
julia> y = [1, 2, 3];
julia> x' * y
* (affine; real)
├─ reshape (affine; real)
│ └─ * (affine; real)
│ ├─ 3×3 SparseArrays.SparseMatrixCSC{Int64, Int64} with 3 stored entries
│ └─ reshape (affine; real)
│ └─ …
└─ [1; 2; 3;;]
+
Base.:+
— MethodBase.:+(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
Base.:+(x::Convex.Value, y::Convex.AbstractExpr)
Base.:+(x::Convex.AbstractExpr, y::Convex.Value)
The addition operator $x + y$.
Examples
Applies to scalar expressions:
julia> x = Variable();
julia> x + 1
+ (affine; real)
├─ real variable (id: 110…477)
└─ [1;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> y = [1, 2, 3];
julia> atom = x + y
+ (affine; real)
├─ 3-element real variable (id: 458…482)
└─ [1; 2; 3;;]
julia> size(atom)
(3, 1)
-
Base.:-
— MethodBase.:-(x::Convex.AbstractExpr)
The univariate negation operator $-x$.
Examples
Applies to scalar expressions:
julia> x = Variable();
julia> -x
Convex.NegateAtom (affine; real)
├─ real variable (id: 161…677)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = -x
Convex.NegateAtom (affine; real)
└─ 3-element real variable (id: 137…541)
julia> size(atom)
(3, 1)
Base.:-
— MethodBase.:-(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
Base.:-(x::Convex.Value, y::Convex.AbstractExpr)
Base.:-(x::Convex.AbstractExpr, y::Convex.Value)
The subtraction operator $x - y$.
Examples
Applies to scalar expressions:
julia> x = Variable();
julia> x - 1
+ (affine; real)
├─ real variable (id: 161…677)
└─ [-1;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> y = [1, 2, 3];
julia> atom = y - x
+ (affine; real)
├─ [1; 2; 3;;]
└─ Convex.NegateAtom (affine; real)
└─ 3-element real variable (id: 242…661)
/
Base.:/
— MethodBase.:/(x::Convex.AbstractExpr, y::Convex.Value)
The binary division operator $\frac{x}{y}$.
Examples
Applies to a scalar expression:
ulia> x = Variable();
julia> x / 2
and element-wise to a matrix:
julia> x = Variable(3);
julia> atom = x / 2
* (affine; real)
├─ 3-element real variable (id: 129…611)
└─ [0.5;;]
julia> size(atom)
(3, 1)
.*
Base.Broadcast.broadcasted
— Methodx::Convex.AbstractExpr .* y::Convex.AbstractExpr
Element-wise multiplication between matrices x
and y
.
Examples
julia> x = Variable(2);
julia> atom = x .* 2
* (affine; real)
├─ 2-element real variable (id: 197…044)
└─ [2;;]
julia> atom = x .* [2, 4]
.* (affine; real)
├─ 2-element real variable (id: 197…044)
└─ [2; 4;;]
julia> size(atom)
(2, 1)
./
Base.Broadcast.broadcasted
— Methodx::Convex.AbstractExpr ./ y::Convex.AbstractExpr
Element-wise division between matrices x
and y
.
Examples
julia> x = Variable(2);
julia> atom = x ./ 2
* (affine; real)
├─ 2-element real variable (id: 875…859)
└─ [0.5;;]
julia> atom = x ./ [2, 4]
.* (affine; real)
├─ 2-element real variable (id: 875…859)
└─ [0.5; 0.25;;]
julia> size(atom)
(2, 1)
.^
Base.Broadcast.broadcasted
— Methodx::Convex.AbstractExpr .^ k::Int
Element-wise exponentiation of x
to the power of k
.
Examples
julia> x = Variable(2);
julia> atom = x .^ 2
qol_elem (convex; positive)
├─ 2-element real variable (id: 131…737)
└─ [1.0; 1.0;;]
julia> size(atom)
(2, 1)
abs
Base.abs
— MethodBase.abs(x::Convex.AbstractExpr)
The epigraph of $|x|$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> abs(x)
abs (convex; positive)
└─ real variable (id: 103…720)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = abs(x)
abs (convex; positive)
└─ 3-element real variable (id: 389…882)
julia> size(atom)
(3, 1)
abs2
Base.abs2
— MethodBase.abs2(x::Convex.AbstractExpr)
The epigraph of $|x|^2$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> abs2(x)
qol_elem (convex; positive)
├─ abs (convex; positive)
│ └─ real variable (id: 319…413)
└─ [1.0;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = abs2(x)
qol_elem (convex; positive)
├─ abs (convex; positive)
│ └─ 3-element real variable (id: 123…996)
└─ [1.0; 1.0; 1.0;;]
julia> size(atom)
(3, 1)
adjoint
Base.adjoint
— MethodLinearAlgebra.adjoint(x::AbstractExpr)
The transpose of the conjugated matrix x
.
Examples
julia> x = ComplexVariable(2, 2);
julia> atom = adjoint(x)
reshape (affine; complex)
└─ * (affine; complex)
├─ 4×4 SparseArrays.SparseMatrixCSC{Int64, Int64} with 4 stored entries
└─ reshape (affine; complex)
└─ conj (affine; complex)
└─ …
julia> size(atom)
(2, 2)
conj
Base.conj
— MethodBase.conj(x::Convex.AbstractExpr)
The complex conjugate of x
.
If x
is real, this function returns x
.
Examples
Applies to a single expression:
julia> x = ComplexVariable();
julia> conj(x)
conj (affine; complex)
└─ complex variable (id: 180…137)
And element-wise to a matrix of expressions:
conj (affine; complex)
└─ complex variable (id: 180…137)
julia> x = ComplexVariable(3);
julia> atom = conj(x)
conj (affine; complex)
└─ 3-element complex variable (id: 104…031)
julia> size(atom)
(3, 1)
conv
Convex.conv
— MethodConvex.conv(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
The convolution between two vectors x
and y
.
Examples
julia> x = Variable(2);
julia> y = [2, 4];
julia> atom = conv(x, y)
* (affine; real)
├─ 3×2 SparseArrays.SparseMatrixCSC{Int64, Int64} with 4 stored entries
└─ 2-element real variable (id: 663…363)
julia> size(atom)
(3, 1)
diag
LinearAlgebra.diag
— FunctionLinearAlgebra.diag(x::Convex.AbstractExpr, k::Int = 0)
Return the k
-th diagonnal of the matrix X
as a column vector.
Examples
Applies to a single square matrix:
julia> x = Variable(2, 2);
julia> atom = diag(x, 0)
diag (affine; real)
└─ 2×2 real variable (id: 724…318)
julia> size(atom)
(2, 1)
julia> atom = diag(x, 1)
diag (affine; real)
└─ 2×2 real variable (id: 147…856)
julia> size(atom)
(1, 1)
diagm
LinearAlgebra.diagm
— MethodLinearAlgebra.diagm(x::Convex.AbstractExpr)
Create a diagonal matrix out of the vector x
.
Examples
julia> x = Variable(2);
julia> atom = diagm(x)
diagm (affine; real)
└─ 2-element real variable (id: 541…968)
julia> size(atom)
(2, 2)
dot
LinearAlgebra.dot
— MethodLinearAlgebra.dot(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
The dot product $x \cdot y$. If x
is complex, it is conjugated.
Examples
julia> x = ComplexVariable(2);
julia> y = [1, 2];
julia> atom = dot(x, y)
sum (affine; complex)
└─ .* (affine; complex)
├─ conj (affine; complex)
│ └─ 2-element complex variable (id: 133…443)
└─ [1; 2;;]
julia> size(atom)
(1, 1)
dotsort
Convex.dotsort
— Methoddotsort(x::Convex.AbstractExpr, y::Convex.Value)
dotsort(x::Convex.Value, y::Convex.AbstractExpr)
Computes dot(sort(x), sort(y))
, where x
or y
is constant.
For example, if x = Variable(6)
and y = [1 1 1 0 0 0]
, this atom computes the sum of the three largest elements of x
.
Examples
julia> x = Variable(4);
julia> atom = dotsort(x, [1, 0, 0, 1])
dotsort (convex; real)
└─ 4-element real variable (id: 128…367)
julia> size(atom)
(1, 1)
eigmax
LinearAlgebra.eigmax
— MethodLinearAlgebra.eigmax(X::Convex.AbstractExpr)
The epigraph of the maximum eigen value of $X$.
Examples
Applies to a single square matrix:
julia> x = Variable(2, 2);
julia> atom = eigmax(x)
eigmin (convex; real)
└─ 2×2 real variable (id: 428…695)
julia> size(atom)
(1, 1)
eigmin
LinearAlgebra.eigmin
— MethodLinearAlgebra.eigmin(X::Convex.AbstractExpr)
The hypograph of the minimum eigen value of $X$.
Examples
Applies to a single square matrix:
julia> x = Variable(2, 2);
julia> atom = eigmin(x)
eigmin (concave; real)
└─ 2×2 real variable (id: 428…695)
julia> size(atom)
(1, 1)
entropy
Convex.entropy
— Methodentropy(x::Convex.AbstractExpr)
The hypograph of $\sum_i -x_i \log x_i$.
Examples
Applies to a matrix of expressions:
julia> x = Variable(3);
julia> atom = entropy(x)
sum (concave; real)
└─ entropy (concave; real)
└─ 3-element real variable (id: 901…778)
julia> size(atom)
(1, 1)
entropy_elementwise
Convex.entropy_elementwise
— Methodentropy_elementwise(x::Convex.AbstractExpr)
The hypograph of $-x \log x$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> entropy_elementwise(x)
entropy (concave; real)
└─ real variable (id: 172…395)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = entropy_elementwise(x)
entropy (concave; real)
└─ 3-element real variable (id: 140…126)
julia> size(atom)
(3, 1)
exp
Base.exp
— MethodBase.exp(x::Convex.AbstractExpr)
The epigraph of $e^x$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> exp(x)
exp (convex; positive)
└─ real variable (id: 103…720)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = exp(x)
exp (convex; positive)
└─ 3-element real variable (id: 389…882)
julia> size(atom)
(3, 1)
geomean
Convex.geomean
— Methodgeomean(x::Convex.AbstractExpr...)
The hypograph of the geometric mean $\sqrt[n]{x_1 \cdot x_2 \cdot \ldots x_n}$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> y = Variable();
julia> geomean(x, y)
geomean (concave; positive)
├─ real variable (id: 163…519)
└─ real variable (id: 107…393)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> y = Variable(3);
julia> atom = geomean(x, y)
geomean (concave; positive)
├─ 3-element real variable (id: 177…782)
└─ 3-element real variable (id: 307…913)
julia> size(atom)
(3, 1)
hcat
Base.hcat
— MethodBase.hcat(args::AbstractExpr...)
Horizontally concatenate args
.
Examples
Applies to a matrix:
julia> x = Variable(2, 2);
julia> atom = hcat(x, x)
hcat (affine; real)
├─ 2×2 real variable (id: 111…376)
└─ 2×2 real variable (id: 111…376)
julia> size(atom)
(2, 4)
You can also use the Julia [x x]
syntax:
julia> x = Variable(2, 2);
julia> atom = [x x]
hcat (affine; real)
├─ 2×2 real variable (id: 111…376)
└─ 2×2 real variable (id: 111…376)
julia> size(atom)
(2, 4)
hinge_loss
Convex.hinge_loss
— Methodhinge_loss(x::Convex.AbstractExpr)
The epigraph of $\max(1 - x, 0)$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> hinge_loss(x)
max (convex; positive)
├─ + (affine; real)
│ ├─ [1;;]
│ └─ Convex.NegateAtom (affine; real)
│ └─ real variable (id: 129…000)
└─ [0;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = hinge_loss(x)
max (convex; positive)
├─ + (affine; real)
│ ├─ * (constant; positive)
│ │ ├─ [1;;]
│ │ └─ [1.0; 1.0; 1.0;;]
│ └─ Convex.NegateAtom (affine; real)
│ └─ 3-element real variable (id: 125…591)
└─ [0;;]
julia> size(atom)
(3, 1)
huber
Convex.huber
— Functionhuber(x::Convex.AbstractExpr, M::Real = 1.0)
The epigraph of the Huber loss function:
\[\begin{cases} x^2 & |x| \le M \\ 2M|x| - M^2 & |x| > M \end{cases}\]
where $M \ge 1$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> huber(x, 2.5)
huber (convex; positive)
└─ real variable (id: 973…369)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = huber(x)
huber (convex; positive)
└─ 3-element real variable (id: 896…728)
julia> size(atom)
(3, 1)
hvcat
Base.hvcat
— MethodBase.hvcat(
rows::Tuple{Vararg{Int}},
args::Union{AbstractExpr,Value}...,
)
Horizontally and vertically concatenate args
in single call.
rows
is the number of arguments to vertically concatenate into each column.
Examples
Applies to a matrix:
To make the matrix:
a b[1] b[2]
c[1] c[2] c[3]
do:
julia> a = Variable();
julia> b = Variable(1, 2);
julia> c = Variable(1, 3);
julia> atom = [a b; c] # Syntactic sugar for: hvcat((2, 1), a, b, c)
vcat (affine; real)
├─ hcat (affine; real)
│ ├─ real variable (id: 429…021)
│ └─ 1×2 real variable (id: 120…326)
└─ hcat (affine; real)
└─ 1×3 real variable (id: 124…615)
julia> size(atom)
(2, 3)
imag
Base.imag
— MethodBase.imag(x::Convex.AbstractExpr)
Return the imaginary component of x
.
Examples
Applies to a single expression:
julia> x = ComplexVariable();
julia> imag(x)
imag (affine; real)
└─ complex variable (id: 407…692)
And element-wise to a matrix of expressions:
julia> x = ComplexVariable(3);
julia> atom = imag(x)
imag (affine; real)
└─ 3-element complex variable (id: 435…057)
julia> size(atom)
(3, 1)
inner_product
Convex.inner_product
— Methodinner_product(x::AbstractExpr, y::AbstractExpr)
The inner product $tr(x^\top y)$ where x
and y
are square matrices.
Examples
julia> x = Variable(2, 2);
julia> y = [1 3; 2 4];
julia> atom = inner_product(x, y)
real (affine; real)
└─ sum (affine; real)
└─ diag (affine; real)
└─ * (affine; real)
├─ …
└─ …
julia> size(atom)
(1, 1)
invpos
Convex.invpos
— Methodinvpos(x::Convex.AbstractExpr)
The epigraph of $\frac{1}{x}$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> invpos(x)
qol_elem (convex; positive)
├─ [1.0;;]
└─ real variable (id: 139…839)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = invpos(x)
qol_elem (convex; positive)
├─ [1.0; 1.0; 1.0;;]
└─ 3-element real variable (id: 133…285)
julia> size(atom)
(3, 1)
kron
Base.kron
— MethodBase.kron(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
The Kronecker (outer) product.
Examples
julia> x = Variable(2);
julia> y = [1 2];
julia> atom = kron(x, y)
vcat (affine; real)
├─ * (affine; real)
│ ├─ index (affine; real)
│ │ └─ 2-element real variable (id: 369…232)
│ └─ [1 2]
└─ * (affine; real)
├─ index (affine; real)
│ └─ 2-element real variable (id: 369…232)
└─ [1 2]
julia> size(atom)
(2, 2)
lieb_ando
Convex.lieb_ando
— Methodlieb_ando(
A::Union{AbstractMatrix,Constant},
B::Union{AbstractMatrix,Constant},
K::Union{AbstractMatrix,Constant},
t::Rational,
)
Returns LinearAlgebra.tr(K' * A^{1-t} * K * B^t)
where A
and B
are positive semidefinite matrices and K
is an arbitrary matrix (possibly rectangular).
lieb_ando(A, B, K, t)
is concave in (A, B)
for t
in [0, 1]
, and convex in (A, B)
for t
in [-1, 0)
or (1, 2]
. K
is a fixed matrix.
Seems numerically unstable when t is on the endpoints of these ranges.
Reference
Ported from CVXQUAD which is based on the paper: "Lieb's concavity theorem, matrix geometric means and semidefinite optimization" by Hamza Fawzi and James Saunderson (arXiv:1512.03401)
Examples
Note that lieb_ando
is implemented as a subproblem, so the returned atom is a Convex.Problem
object. The Problem
atom can still be used as a regular 1x1
atom in other expressions.
julia> A = Semidefinite(2, 2);
julia> B = Semidefinite(3, 3);
julia> K = [1 2 3; 4 5 6];
julia> atom = lieb_ando(A, B, K, 1 // 2)
Problem statistics
problem is DCP : true
number of variables : 3 (49 scalar elements)
number of constraints : 4 (157 scalar elements)
number of coefficients : 76
number of atoms : 26
Solution summary
termination status : OPTIMIZE_NOT_CALLED
primal status : NO_SOLUTION
dual status : NO_SOLUTION
Expression graph
maximize
└─ real (affine; real)
└─ sum (affine; real)
└─ diag (affine; real)
└─ …
subject to
├─ GeometricMeanHypoConeSquare constraint (convex)
│ └─ vcat (affine; real)
│ ├─ reshape (affine; real)
│ │ └─ …
│ ├─ reshape (affine; real)
│ │ └─ …
│ └─ reshape (affine; real)
│ └─ …
├─ PSD constraint (convex)
│ └─ 6×6 real variable (id: 173…902)
├─ PSD constraint (convex)
│ └─ 6×6 real variable (id: 173…902)
⋮
julia> size(atom)
(1, 1)
log
Base.log
— MethodBase.log(x::Convex.AbstractExpr)
The hypograph of $\log(x)$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> log(x)
log (concave; real)
└─ real variable (id: 103…720)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = log(x)
log (concave; real)
└─ 3-element real variable (id: 161…499)
julia> size(atom)
(3, 1)
log_perspective
Convex.log_perspective
— Methodlog_perspective(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
The hypograph the perspective of of the log function: $\sum y_i*\log \frac{x_i}{y_i}$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> y = Variable();
julia> log_perspective(x, y)
Convex.NegateAtom (concave; real)
└─ relative_entropy (convex; real)
├─ real variable (id: 136…971)
└─ real variable (id: 131…344)
And to a matrix of expressions:
julia> x = Variable(3);
julia> y = Variable(3);
julia> atom = log_perspective(x, y)
Convex.NegateAtom (concave; real)
└─ relative_entropy (convex; real)
├─ 3-element real variable (id: 854…248)
└─ 3-element real variable (id: 111…174)
julia> size(atom)
(1, 1)
logdet
LinearAlgebra.logdet
— MethodLinearAlgebra.logdet(X::Convex.AbstractExpr)
The hypograph of $\log(\det(X))$.
Examples
Applies to a single matrix expression:
julia> X = Variable(2, 2);
julia> atom = logdet(X)
logdet (concave; real)
└─ 2×2 real variable (id: 159…883)
julia> size(atom)
(1, 1)
logisticloss
Convex.logisticloss
— Methodlogisticloss(x::Convex.AbstractExpr)
Reformulation for epigraph of the logistic loss: $\sum_i \log(e^x_i + 1)$.
This reformulation uses logsumexp
.
Examples
Applies to a single expression:
julia> x = Variable();
julia> logisticloss(x)
logsumexp (convex; real)
└─ vcat (affine; real)
├─ real variable (id: 444…892)
└─ [0;;]
And to a matrix of expressions:
julia> x = Variable(3);
julia> atom = logisticloss(x)
+ (convex; real)
├─ logsumexp (convex; real)
│ └─ vcat (affine; real)
│ ├─ index (affine; real)
│ │ └─ …
│ └─ [0;;]
├─ logsumexp (convex; real)
│ └─ vcat (affine; real)
│ ├─ index (affine; real)
│ │ └─ …
│ └─ [0;;]
└─ logsumexp (convex; real)
└─ vcat (affine; real)
├─ index (affine; real)
│ └─ …
└─ [0;;]
julia> size(atom)
(1, 1)
logsumexp
Convex.logsumexp
— Methodlogsumexp(x::Convex.AbstractExpr)
The epigraph of $\log\left(\sum_i e^{x_i}\right)$.
Examples
Applies to a single expression:
julia> x = Variable(2, 3);
julia> atom = logsumexp(x)
logsumexp (convex; real)
└─ 2×3 real variable (id: 121…604)
julia> size(atom)
(1, 1)
julia> atom = logsumexp(x; dims = 1)
logsumexp (convex; real)
└─ 2×3 real variable (id: 121…604)
julia> size(atom)
(1, 3)
julia> atom = logsumexp(x; dims = 2)
logsumexp (convex; real)
└─ 2×3 real variable (id: 121…604)
julia> size(atom)
(2, 1)
matrixfrac
Convex.matrixfrac
— Methodmatrixfrac(x::AbstractExpr, P::AbstractExpr)
The epigraph of $x^\top P^{-1} x$.
Examples
julia> x = Variable(2);
julia> P = Variable(2, 2);
julia> atom = matrixfrac(x, P)
matrixfrac (convex; positive)
├─ 2-element real variable (id: 139…388)
└─ 2×2 real variable (id: 126…414)
julia> size(atom)
(1, 1)
max
Base.max
— MethodBase.max(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
Base.max(x::Convex.AbstractExpr, y::Convex.Value)
Base.max(x::Convex.Value, y::Convex.AbstractExpr)
The hypograph of $max(x, y)$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> max(x, 1)
max (convex; real)
├─ real variable (id: 183…974)
└─ [1;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> y = [1, 2, 3];
julia> atom = max(x, y)
max (convex; real)
├─ 3-element real variable (id: 153…965)
└─ [1; 2; 3;;]
julia> size(atom)
(3, 1)
maximum
Base.maximum
— MethodBase.maximum(x::Convex.AbstractExpr)
The hypograph of $max(x...)$.
Examples
Applies to a matrix expression:
julia> x = Variable(3);
julia> atom = maximum(x)
maximum (convex; real)
└─ 3-element real variable (id: 159…219)
julia> size(atom)
(1, 1)
min
Base.min
— MethodBase.min(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
Base.min(x::Convex.Value, y::Convex.AbstractExpr)
Base.min(x::Convex.AbstractExpr, y::Convex.Value)
The epigraph of $min(x, y)$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> min(x, 1)
min (concave; real)
├─ real variable (id: 183…974)
└─ [1;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> y = [1, 2, 3];
julia> atom = min(x, y)
min (concave; real)
├─ 3-element real variable (id: 153…965)
└─ [1; 2; 3;;]
julia> size(atom)
(3, 1)
minimum
Base.minimum
— MethodBase.minimum(x::Convex.AbstractExpr)
The epigraph of $min(x...)$.
Examples
Applies to a matrix expression:
julia> x = Variable(3);
julia> atom = minimum(x)
minimum (convex; real)
└─ 3-element real variable (id: 159…219)
julia> size(atom)
(1, 1)
neg
Convex.neg
— Methodneg(x::Convex.AbstractExpr)
The epigraph of $\max(-x, 0)$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> neg(x)
max (convex; positive)
├─ Convex.NegateAtom (affine; real)
│ └─ real variable (id: 467…111)
└─ [0;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = neg(x)
max (convex; positive)
├─ Convex.NegateAtom (affine; real)
│ └─ 3-element real variable (id: 224…439)
└─ [0;;]
julia> size(atom)
(3, 1)
norm
LinearAlgebra.norm
— Functionnorm(x::AbstractExpr, p::Real = 2)
Computes the p
-norm ‖x‖ₚ = (∑ᵢ |xᵢ|^p)^(1/p)
of a vector expression x
.
Matrices are vectorized (i.e., norm(x)
is the same as norm(vec(x))
.)
The return value depends on the value of p
. Specialized cases are used for p = 1
, p = 2
, and p = Inf
.
Examples
julia> x = Variable(2);
julia> atom = norm(x, 1)
sum (convex; positive)
└─ abs (convex; positive)
└─ 2-element real variable (id: 779…899)
julia> size(atom)
(1, 1)
julia> norm(x, 2)
norm2 (convex; positive)
└─ 2-element real variable (id: 779…899)
julia> norm(x, Inf)
maximum (convex; positive)
└─ abs (convex; positive)
└─ 2-element real variable (id: 779…899)
julia> norm(x, 3 // 2)
rationalnorm (convex; positive)
└─ 2-element real variable (id: 779…899)
norm2
LinearAlgebra.norm2
— MethodLinearAlgebra.norm2(x::Convex.AbstractExpr)
The epigraph of the 2-norm $||x||_2$.
Examples
Applies to a matrix of expressions:
julia> x = Variable(3);
julia> atom = norm2(x)
norm2 (convex; positive)
└─ 3-element real variable (id: 162…975)
julia> size(atom)
(3, 1)
And to a complex:
julia> y = ComplexVariable(3);
julia> atom = norm2(y)
norm2 (convex; positive)
└─ vcat (affine; real)
├─ real (affine; real)
│ └─ 3-element complex variable (id: 120…942)
└─ imag (affine; real)
└─ 3-element complex variable (id: 120…942)
julia> size(atom)
(1, 1)
nuclearnorm
Convex.nuclearnorm
— Methodnuclearnorm(x::Convex.AbstractExpr)
The epigraph of the nuclear norm $||X||_*$, which is the sum of the singular values of $X$.
Examples
Applies to a real-valued matrix:
julia> x = Variable(2, 2);
julia> atom = nuclearnorm(x)
nuclearnorm (convex; positive)
└─ 2×2 real variable (id: 106…758)
julia> size(atom)
(1, 1)
julia> y = ComplexVariable(2, 2);
julia> atom = nuclearnorm(y)
nuclearnorm (convex; positive)
└─ 2×2 complex variable (id: 577…313)
julia> size(atom)
(1, 1)
opnorm
LinearAlgebra.opnorm
— FunctionLinearAlgebra.opnorm(x::Convex.AbstractExpr, p::Real = 2)
The epigraph of the matrix norm $||X||_p$.
Examples
Applies to a real- or complex-valued matrix:
julia> x = Variable(2, 2);
julia> atom = LinearAlgebra.opnorm(x, 1)
maximum (convex; positive)
└─ * (convex; positive)
├─ [1.0 1.0]
└─ abs (convex; positive)
└─ 2×2 real variable (id: 106…758)
julia> atom = LinearAlgebra.opnorm(x, 2)
opnorm (convex; positive)
└─ 2×2 real variable (id: 106…758)
julia> atom = LinearAlgebra.opnorm(x, Inf)
maximum (convex; positive)
└─ * (convex; positive)
├─ abs (convex; positive)
│ └─ 2×2 real variable (id: 106…758)
└─ [1.0; 1.0;;]
julia> y = ComplexVariable(2, 2);
julia> atom = maximum (convex; positive)
└─ * (convex; positive)
├─ abs (convex; positive)
│ └─ 2×2 complex variable (id: 116…943)
└─ [1.0; 1.0;;]
julia> size(atom)
(1, 1)
partialtrace
Convex.partialtrace
— Methodpartialtrace(x, sys::Int, dims::Vector)
Returns the partial trace of x
over the sys
th system, where dims
is a vector of integers encoding the dimensions of each subsystem.
partialtranspose
Convex.partialtranspose
— Methodpartialtranspose(x, sys::Int, dims::Vector)
Returns the partial transpose of x
over the sys
th system, where dims
is a vector of integers encoding the dimensions of each subsystem.
pos
Convex.pos
— Methodpos(x::Convex.AbstractExpr)
The epigraph of $\max(x, 0)$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> pos(x)
max (convex; positive)
├─ real variable (id: 467…111)
└─ [0;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = pos(x)
max (convex; positive)
├─ 3-element real variable (id: 154…809)
└─ [0;;]
julia> size(atom)
(3, 1)
qol_elementwise
Convex.qol_elementwise
— Methodqol_elementwise(x::AbstractExpr, y::AbstractExpr)
The elementwise epigraph of $\frac{x^2}{y}$.
Examples
julia> x = Variable(3);
julia> y = Variable(3, Positive());
julia> atom = qol_elementwise(x, y)
qol_elem (convex; positive)
├─ 3-element real variable (id: 155…648)
└─ 3-element positive variable (id: 227…080)
julia> size(atom)
(3, 1)
quadform
Convex.quadform
— Methodquadform(x::AbstractExpr, A::AbstractExpr; assume_psd=false)
Represents $x^\top A x$ where either:
x
is a vector-valued variable andA
is a positive semidefinite or negative semidefinite matrix (and in particular Hermitian or real symmetric). Ifassume_psd=true
, thenA
will be assumed to be positive semidefinite. Otherwise,Convex._is_psd
will be used to check ifA
is positive semidefinite or negative semidefinite.- or
A
is a matrix-valued variable andx
is a vector.
Examples
julia> x = Variable(2);
julia> A = [1 0; 0 1]
2×2 Matrix{Int64}:
1 0
0 1
julia> atom = quadform(x, A)
* (convex; positive)
├─ [1;;]
└─ qol_elem (convex; positive)
├─ norm2 (convex; positive)
│ └─ * (affine; real)
│ ├─ …
│ └─ …
└─ [1.0;;]
julia> size(atom)
(1, 1)
julia> x = [1, 2]
julia> A = Variable(2, 2);
julia> atom = quadform(x, A)
* (affine; real)
├─ * (affine; real)
│ ├─ [1 2]
│ └─ 2×2 real variable (id: 111…794)
└─ [1; 2;;]
julia> size(atom)
(1, 1)
quadoverlin
Convex.quadoverlin
— Methodquadoverlin(x::AbstractExpr, y::AbstractExpr)
The epigraph of $\frac{||x||_2^2}{y}$.
Examples
julia> x = Variable(3);
julia> y = Variable(Positive());
julia> atom = quadoverlin(x, y)
qol (convex; positive)
├─ 3-element real variable (id: 868…883)
└─ positive variable (id: 991…712)
julia> size(atom)
(1, 1)
quantum_entropy
Convex.quantum_entropy
— Functionquantum_entropy(X::AbstractExpr, m::Integer, k::Integer)
quantum_entropy
returns -LinearAlgebra.tr(X*log(X))
where X
is a positive semidefinite.
Note this function uses logarithm base e, not base 2, so return value is in units of nats, not bits.
Quantum entropy is concave. This function implements the semidefinite programming approximation given in the reference below. Parameters m
and k
control the accuracy of this approximation: m
is the number of quadrature nodes to use and k
the number of square-roots to take. See reference for more details.
The implementation uses the expression
\[H(X) = -tr(D_{op}(X||I))\]
where $D_{op}$ is the operator relative entropy:
\[D_{op}(X||Y) = X^{1/2}*logm(X^{1/2} Y^{-1} X^{1/2})*X^{1/2}\]
Reference
Ported from CVXQUAD which is based on the paper: "Lieb's concavity theorem, matrix geometric means and semidefinite optimization" by Hamza Fawzi and James Saunderson (arXiv:1512.03401)
Examples
Applies to a matrix:
julia> X = Variable(2, 2);
julia> atom = quantum_entropy(X)
quantum_entropy (concave; positive)
└─ 2×2 real variable (id: 700…694)
julia> size(atom)
(1, 1)
quantum_relative_entropy
Convex.quantum_relative_entropy
— Functionquantum_relative_entropy(
A::AbstractExpr,
B::AbstractExpr,
m::Integer,
k::Integer,
)
quantum_relative_entropy
returns LinearAlgebra.tr(A*(log(A)-log(B)))
where A
and B
are positive semidefinite matrices.
Note this function uses logarithm base e, not base 2, so return value is in units of nats, not bits.
Quantum relative entropy is convex (jointly) in (A, B). This function implements the semidefinite programming approximation given in the reference below. Parameters m and k control the accuracy of this approximation: m
is the number of quadrature nodes to use and k
the number of square-roots to take. See reference for more details.
Implementation uses the expression
\[D(A||B) = e'*D_{op} (A \otimes I || I \otimes B) )*e\]
where $D_{op}$ is the operator relative entropy and e = vec(Matrix(I, n, n))
.
Reference
Ported from CVXQUAD which is based on the paper: "Lieb's concavity theorem, matrix geometric means and semidefinite optimization" by Hamza Fawzi and James Saunderson (arXiv:1512.03401)
Examples
julia> A = Variable(2, 2);
julia> B = Variable(2, 2);
julia> atom = quantum_relative_entropy(A, B)
quantum_relative_entropy (convex; positive)
├─ 2×2 real variable (id: 144…849)
└─ 2×2 real variable (id: 969…693)
julia> size(atom)
(1, 1)
rationalnorm
Convex.rationalnorm
— Methodrationalnorm(x::AbstractExpr, k::Rational{Int})
The epigraph of ||x||_k
.
Examples
Applies to a single matrix:
julia> x = Variable(2);
julia> atom = rationalnorm(x, 3 // 2)
rationalnorm (convex; positive)
└─ 2-element real variable (id: 182…293)
julia> size(atom)
(1, 1)
real
Base.real
— MethodBase.real(x::Convex.AbstractExpr)
Return the real component of x
.
Examples
Applies to a single expression:
julia> x = ComplexVariable();
julia> real(x)
real (affine; real)
└─ complex variable (id: 407…692)
And element-wise to a matrix of expressions:
julia> x = ComplexVariable(3);
julia> atom = real(x)
real (affine; real)
└─ 3-element complex variable (id: 435…057)
julia> size(atom)
(3, 1)
relative_entropy
Convex.relative_entropy
— Methodrelative_entropy(x::Convex.AbstractExpr, y::Convex.AbstractExpr)
The epigraph of $\sum x_i*\log \frac{x_i}{y_i}$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> y = Variable();
julia> relative_entropy(x, y)
relative_entropy (convex; real)
├─ real variable (id: 124…372)
└─ real variable (id: 409…346)
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> y = Variable(3);
julia> atom = relative_entropy(x, y)
relative_entropy (convex; real)
├─ 3-element real variable (id: 906…671)
└─ 3-element real variable (id: 118…912)
julia> size(atom)
(1, 1)
reshape
Base.reshape
— MethodBase.reshape(x::AbstractExpr, m::Int, n::Int)
Reshapes the expression x
into a matrix with m
rows and n
columns.
Examples
Applies to a matrix:
julia> x = Variable(6, 1);
julia> size(x)
(6, 1)
julia> atom = reshape(x, 2, 3)
reshape (affine; real)
└─ 6-element real variable (id: 103…813)
julia> size(atom)
(2, 3)
rootdet
Convex.rootdet
— MethodConvex.rootdet(X::Convex.AbstractExpr)
The hypograph of $\det(X)^{\frac{1}{n}}$, where $n$ is the side-dimension of the square matrix $X$.
Examples
Applies to a single matrix expression:
julia> X = Variable(2, 2);
julia> atom = rootdet(X)
rootdet (concave; real)
└─ 2×2 real variable (id: 159…883)
julia> size(atom)
(1, 1)
sigmamax
Convex.sigmamax
— Methodsigmamax(x::Convex.AbstractExpr)
The epigraph of the spectral norm $||X||_2$, which is the maximum of the singular values of $X$.
Examples
Applies to a real- or complex-valued matrix:
julia> x = Variable(2, 2);
julia> atom = sigmamax(x)
opnorm (convex; positive)
└─ 2×2 real variable (id: 106…758)
julia> size(atom)
(1, 1)
julia> y = ComplexVariable(2, 2);
julia> atom = sigmamax(y)
opnorm (convex; positive)
└─ 2×2 complex variable (id: 577…313)
julia> size(atom)
(1, 1)
sqrt
Base.sqrt
— MethodBase.sqrt(x::Convex.AbstractExpr)
The hypograph of $\sqrt x$.
Examples
Applies to a single expression:
julia> x = Variable();
julia> sqrt(x)
geomean (concave; positive)
├─ real variable (id: 576…546)
└─ [1.0;;]
And element-wise to a matrix of expressions:
julia> x = Variable(3);
julia> atom = sqrt(x)
geomean (concave; positive)
├─ 3-element real variable (id: 181…583)
└─ [1.0; 1.0; 1.0;;]
julia> size(atom)
(3, 1)
square
Convex.square
— Methodsquare(x::AbstractExpr)
The epigraph of $x^2$.
Examples
Applies elementwise to a matrix
julia> x = Variable(3);
julia> atom = square(x)
qol_elem (convex; positive)
├─ 3-element real variable (id: 438…681)
└─ [1.0; 1.0; 1.0;;]
julia> size(atom)
(3, 1)
sum
Base.sum
— MethodBase.sum(x::Convex.AbstractExpr; dims = :)
Sum x
, optionally along a dimension dims
.
Examples
Sum all elements in an expression:
julia> x = Variable(2, 2);
julia> atom = sum(x)
sum (affine; real)
└─ 2×2 real variable (id: 263…449)
julia> size(atom)
(1, 1)
Sum along the first dimension, creating a row vector:
julia> x = Variable(2, 2);
julia> atom = sum(x; dims = 1)
* (affine; real)
├─ [1.0 1.0]
└─ 2×2 real variable (id: 143…826)
julia> size(atom)
(1, 2)
Sum along the second dimension, creating a columnn vector:
julia> atom = sum(x; dims = 2)
* (affine; real)
├─ 2×2 real variable (id: 143…826)
└─ [1.0; 1.0;;]
julia> size(atom)
(2, 1)
sumlargest
Convex.sumlargest
— Methodsumlargest(x::Convex.AbstractExpr, k::Int)
Sum the k
largest values of x
.
Examples
Applies to a matrix:
julia> x = Variable(3, 3);
julia> atom = sumlargest(x, 2)
sumlargest (convex; real)
├─ 3×3 real variable (id: 833…482)
└─ [2;;]
julia> size(atom)
(1, 1)
sumlargesteigs
Convex.sumlargesteigs
— Methodsumlargesteigs(x::Convex.AbstractExpr, k::Int)
Sum the k
largest eigen values of x
.
Examples
Applies to a matrix:
julia> x = Variable(3, 3);
julia> atom = sumlargesteigs(x, 2)
sumlargesteigs (convex; real)
├─ 3×3 real variable (id: 833…482)
└─ [2;;]
julia> size(atom)
(1, 1)
sumsmallest
Convex.sumsmallest
— Methodsumsmallest(x::Convex.AbstractExpr, k::Int)
Sum the k
smallest values of x
.
Examples
Applies to a matrix:
julia> x = Variable(3, 3);
julia> atom = sumsmallest(x, 2)
Convex.NegateAtom (concave; real)
└─ sumlargest (convex; real)
└─ Convex.NegateAtom (affine; real)
└─ 3×3 real variable (id: 723…082)
julia> size(atom)
(1, 1)
sumsquares
Convex.sumsquares
— Methodsumsquares(x::AbstractExpr)
The epigraph of $||x||_2^2$.
Examples
Applies to a single matrix
julia> x = Variable(3);
julia> atom = sumsquares(x)
qol (convex; positive)
├─ 3-element real variable (id: 125…181)
└─ [1;;]
julia> size(atom)
(1, 1)
tr
LinearAlgebra.tr
— MethodLinearAlgebra.tr(x::AbstractExpr)
The trace of the matrix x
.
Examples
julia> x = Variable(2, 2);
julia> atom = tr(x)
sum (affine; real)
└─ diag (affine; real)
└─ 2×2 real variable (id: 844…180)
julia> size(atom)
(1, 1)
trace_logm
Convex.trace_logm
— Functiontrace_logm(
X::Convex.AbstractExpr,
C::AbstractMatrix,
m::Integer = 3,
k::Integer = 3,
)
trace_logm(X, C)
returns LinearAlgebra.tr(C*logm(X))
where X
and C
are positive definite matrices and C
is constant.
trace_logm
is concave in X
.
This function implements the semidefinite programming approximation given in the reference below. Parameters m
and k
control the accuracy of the approximation: m
is the number of quadrature nodes to use and k
is the number of square-roots to take. See reference for more details.
Implementation uses the expression
\[tr(C \times logm(X)) = -tr(C \times D_{op}(I||X))\]
where D_{op}
is the operator relative entropy:
\[D_{op}(X||Y) = X^{1/2}*logm(X^{1/2} Y^{-1} X^{1/2})*X^{1/2}\]
Reference
Ported from CVXQUAD which is based on the paper: "Lieb's concavity theorem, matrix geometric means and semidefinite optimization" by Hamza Fawzi and James Saunderson (arXiv:1512.03401)
Examples
Applies to a matrix:
julia> X = Variable(2, 2);
julia> C = [1 0; 0 1];
julia> atom = trace_logm(X, C)
trace_logm (concave; real)
└─ 2×2 real variable (id: 608…362)
julia> size(atom)
(1, 1)
trace_mpower
Convex.trace_mpower
— Methodtrace_mpower(A::Convex.AbstractExpr, t::Rational, C::AbstractMatrix)
trace_mpower(A, t, C)
returns LinearAlgebra.tr(C*A^t)
where A
and C
are positive definite matrices, C
is constant and t
is a rational in [-1, 2]
.
When t
is in [0, 1]
, trace_mpower(A, t, C)
is concave in A
(for fixed positive semidefinite matrix C
) and convex for t
in [-1, 0)
or (1, 2]
.
Reference
Ported from CVXQUAD which is based on the paper: "Lieb's concavity theorem, matrix geometric means and semidefinite optimization" by Hamza Fawzi and James Saunderson (arXiv:1512.03401)
Examples
Applies to a matrix:
julia> A = Variable(2, 2);
julia> C = [1 0; 0 1];
julia> atom = trace_mpower(A, 1 // 2, C)
trace_mpower (concave; real)
└─ 2×2 real variable (id: 150…626)
julia> size(atom)
(1, 1)
transpose
Base.transpose
— MethodLinearAlgebra.transpose(x::AbstractExpr)
The transpose of the matrix x
.
Examples
julia> x = Variable(2, 2);
julia> atom = transpose(x)
reshape (affine; real)
└─ * (affine; real)
├─ 4×4 SparseArrays.SparseMatrixCSC{Int64, Int64} with 4 stored entries
└─ reshape (affine; real)
└─ 2×2 real variable (id: 151…193)
julia> size(atom)
(2, 2)
vcat
Base.vcat
— MethodBase.vcat(args::AbstractExpr...)
Vertically concatenate args
.
Examples
Applies to a matrix:
julia> x = Variable(2, 2);
julia> atom = vcat(x, x)
vcat (affine; real)
├─ 2×2 real variable (id: 111…376)
└─ 2×2 real variable (id: 111…376)
julia> size(atom)
(4, 2)
You can also use the Julia [x; x]
syntax:
julia> x = Variable(2, 2);
julia> atom = [x; x]
vcat (affine; real)
├─ 2×2 real variable (id: 111…376)
└─ 2×2 real variable (id: 111…376)
julia> size(atom)
(4, 2)
vec
Base.vec
— MethodBase.vec(x::AbstractExpr)
Reshapes the expression x
into a column vector.
Examples
Applies to a matrix:
julia> x = Variable(2, 2);
julia> atom = vec(x)
reshape (affine; real)
└─ 2×2 real variable (id: 115…295)
julia> size(atom)
(4, 1)