Cones module

See the Cone interface and predefined cones section.

Array utilities

Hypatia.Cones.vec_copyto!Function
vec_copyto!(v1, v2)

Copy a vector in-place.

source
vec_copyto!(rvec, cvec)

Copy a complex vector to a real vector in-place.

source
vec_copyto!(cvec, rvec)

Copy a real vector to a complex vector in-place.

source
Hypatia.Cones.svec_lengthFunction
svec_length(side)

The dimension of the vectorized triangle of a real symmetric matrix with side dimension side::Int.

source
svec_length(_, side)

The dimension of the real vectorized triangle of a real symmetric or complex Hermitian matrix with side dimension side::Int.

source
Hypatia.Cones.svec_sideFunction
svec_side(len)

The side dimension of a real symmetric matrix with vectorized triangle length len::Int.

source
svec_side(_, len)

The side dimension of a real symmetric or complex Hermitian matrix with real vectorized triangle length len::Int.

source
Hypatia.Cones.svec_idxFunction
svec_idx(row, col)

The index in the vectorized triangle of a symmetric matrix for element (row::Int, col::Int).

source
Hypatia.Cones.block_idxsFunction
block_idxs(incr, block)

The indices corresponding to block block::Int in a vector of blocks with equal length incr::Int.

source
Hypatia.Cones.scale_svec!Function
scale_svec!(arr, scal)

Rescale the elements corresponding to off-diagonals in arr::AbstractVecOrMat, with scaling scal::Real, for real symmetric matrices.

source
Hypatia.Cones.scale_svec_incr!Function
scale_svec_incr!(arr, scal, incr)

Rescale the elements corresponding to off-diagonals in arr::AbstractVecOrMat, with scaling scal::Real and block increment incr::Int, for real symmetric matrices.

source
Hypatia.Cones.smat_to_svec!Function
smat_to_svec!(vec, mat, rt2)

Copy a real symmetric matrix upper triangle to a svec-scaled vector in-place.

source
smat_to_svec!(vec, mat, rt2)

Copy a complex Hermitian matrix upper triangle to a svec-scaled real vector in-place.

source
Hypatia.Cones.svec_to_smat!Function
svec_to_smat!(mat, vec, rt2)

Copy a svec-scaled vector to a real symmetric matrix upper triangle in-place.

source
svec_to_smat!(mat, vec, rt2)

Copy a svec-scaled real vector to a complex Hermitian matrix upper triangle in-place.

source

Cone oracles

Hypatia.Cones.is_feasFunction
is_feas(cone)

Returns true if and only if the currently-loaded primal point is strictly feasible for the cone.

source
Hypatia.Cones.is_dual_feasFunction
is_dual_feas(cone)

Returns false only if the currently-loaded dual point is outside the interior of the cone's dual cone.

source
Hypatia.Cones.gradFunction
grad(cone)

The gradient of the cone's barrier function at the currently-loaded primal point.

source
Hypatia.Cones.hessFunction
hess(cone)

The Hessian (symmetric positive definite) of the cone's barrier function at the currently-loaded primal point.

source
Hypatia.Cones.inv_hessFunction
inv_hess(cone)

The inverse Hessian (symmetric positive definite) of the cone's barrier function at the currently-loaded primal point.

source
Hypatia.Cones.hess_prod!Function
hess_prod!(prod, arr, cone)

Compute the product of the Hessian of the cone's barrier function at the currently-loaded primal point with a vector or array, in-place.

source
Hypatia.Cones.inv_hess_prod!Function
inv_hess_prod!(prod, arr, cone)

Compute the product of the inverse Hessian of the cone's barrier function at the currently-loaded primal point with a vector or array, in-place.

source
Hypatia.Cones.dder3Function
dder3(cone, dir)

Compute the third-order directional derivative, in the direction dir, the cone's barrier function at the currently-loaded primal point.

source

Predefined cone types

Hypatia.Cones.NonnegativeType
mutable struct Nonnegative{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Nonnegative cone of dimension dim.

Nonnegative{T}(dim::Int)
source
Hypatia.Cones.PosSemidefTriType
mutable struct PosSemidefTri{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Real symmetric or complex Hermitian positive semidefinite cone of dimension dim in svec (scaled upper triangle) format.

PosSemidefTri{T, R}(dim::Int)
source
Hypatia.Cones.DoublyNonnegativeTriType
mutable struct DoublyNonnegativeTri{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Real symmetric doubly nonnegative cone (intersection of nonnegative and positive semidefinite cones) of dimension dim in svec format.

DoublyNonnegativeTri{T}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.PosSemidefTriSparseType
mutable struct PosSemidefTriSparse{I<:Hypatia.Cones.PSDSparseImpl, T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Real symmetric or complex Hermitian sparse positive semidefinite cone of side dimension side and sparse lower triangle row and column indices rows, cols in svec format. Note all diagonal elements must be present.

PosSemidefTriSparse{T, R}(side::Int, rows::Vector{Int}, cols::Vector{Int}, use_dual::Bool = false)
source
Hypatia.Cones.LinMatrixIneqType
mutable struct LinMatrixIneq{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Linear matrix inequality cone parametrized by list of real symmetric or complex Hermitian matrices mats of equal dimension.

LinMatrixIneq{T}(mats::Vector, use_dual::Bool = false)
source
Hypatia.Cones.EpiNormInfType
mutable struct EpiNormInf{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of real or complex infinity norm cone of dimension dim.

EpiNormInf{T, R}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.EpiNormEuclType
mutable struct EpiNormEucl{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of Euclidean norm (AKA second-order) cone of dimension dim.

EpiNormEucl{T}(dim::Int)
source
Hypatia.Cones.EpiPerSquareType
mutable struct EpiPerSquare{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of perspective function of halved squared Euclidean norm (AKA rotated second-order) cone of dimension dim.

EpiPerSquare{T}(dim::Int)
source
Hypatia.Cones.EpiNormSpectralTriType
mutable struct EpiNormSpectralTri{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of real symmetric or complex Hermitian matrix spectral norm (i.e. maximum absolute value of eigenvalues) cone of dimension dim in svec format.

EpiNormSpectralTri{T, R}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.EpiNormSpectralType
mutable struct EpiNormSpectral{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of real or complex matrix spectral norm (i.e. maximum singular value) for a matrix (stacked column-wise) of nrows rows and ncols columns with nrows ≤ ncols.

EpiNormSpectral{T, R}(nrows::Int, ncols::Int, use_dual::Bool = false)
source
Hypatia.Cones.MatrixEpiPerSquareType
mutable struct MatrixEpiPerSquare{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Matrix epigraph of perspective function of real or complex matrix outer product for a matrix (stacked column-wise) of nrows rows and ncols columns with nrows ≤ ncols.

MatrixEpiPerSquare{T, R}(nrows::Int, ncols::Int, use_dual::Bool = false)
source
Hypatia.Cones.GeneralizedPowerType
mutable struct GeneralizedPower{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Generalized power cone parametrized by powers α in the unit simplex and dimension d of the normed variables.

GeneralizedPower{T}(α::Vector{T}, d::Int, use_dual::Bool = false)
source
Hypatia.Cones.HypoPowerMeanType
mutable struct HypoPowerMean{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Hypograph of weighted power mean cone parametrized by powers α in the unit simplex.

HypoPowerMean{T}(α::Vector{T}, use_dual::Bool = false)
source
Hypatia.Cones.HypoGeoMeanType
mutable struct HypoGeoMean{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Hypograph of geometric mean cone of dimension dim.

HypoGeoMean{T}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.HypoRootdetTriType
mutable struct HypoRootdetTri{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Hypograph of real symmetric or complex Hermitian root-determinant cone of dimension dim in svec format.

HypoRootdetTri{T, R}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.HypoPerLogType
mutable struct HypoPerLog{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Hypograph of perspective function of sum-log cone of dimension dim.

HypoPerLog{T}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.HypoPerLogdetTriType
mutable struct HypoPerLogdetTri{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Hypograph of perspective function of real symmetric or complex Hermitian log-determinant cone of dimension dim in svec format.

HypoPerLogdetTri{T, R}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.EpiPerSepSpectralType
mutable struct EpiPerSepSpectral{Q<:Hypatia.Cones.ConeOfSquares, T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of perspective function of a convex separable spectral function h over a cone of squares Q on a Jordan algebra with rank d.

EpiPerSepSpectral{Q, T}(h::Hypatia.Cones.SepSpectralFun, d::Int, use_dual::Bool = false)
source
Hypatia.Cones.EpiRelEntropyType
mutable struct EpiRelEntropy{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of vector relative entropy cone of dimension dim.

EpiRelEntropy{T}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.EpiTrRelEntropyTriType
mutable struct EpiTrRelEntropyTri{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Epigraph of matrix relative entropy cone of dimension dim in svec format.

EpiTrRelEntropyTri{T}(dim::Int, use_dual::Bool = false)
source
Hypatia.Cones.WSOSInterpNonnegativeType
mutable struct WSOSInterpNonnegative{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.Cone{T<:Real}

Interpolant-basis weighted sum-of-squares polynomial cone of dimension U, for real or real-valued complex polynomials, parametrized by vector of matrices Ps derived from interpolant basis and polynomial domain constraints.

WSOSInterpNonnegative{T, R}(U::Int, Ps::Vector{Matrix{R}}, use_dual::Bool = false)
source
Hypatia.Cones.WSOSInterpPosSemidefTriType
mutable struct WSOSInterpPosSemidefTri{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Interpolant-basis weighted sum-of-squares polynomial (of dimension U) positive semidefinite matrix (of side dimension R) cone, parametrized by vector of matrices Ps derived from interpolant basis and polynomial domain constraints.

WSOSInterpPosSemidefTri{T}(R::Int, U::Int, Ps::Vector{Matrix{T}}, use_dual::Bool = false)
source
Hypatia.Cones.WSOSInterpEpiNormEuclType
mutable struct WSOSInterpEpiNormEucl{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Interpolant-basis weighted sum-of-squares polynomial (of dimension U) epigraph of ℓ₂ norm (of dimension R) cone, parametrized by vector of matrices Ps derived from interpolant basis and polynomial domain constraints.

WSOSInterpEpiNormEucl{T}(R::Int, U::Int, Ps::Vector{Matrix{T}}, use_dual::Bool = false)
source
Hypatia.Cones.WSOSInterpEpiNormOneType
mutable struct WSOSInterpEpiNormOne{T<:Real} <: Hypatia.Cones.Cone{T<:Real}

Interpolant-basis weighted sum-of-squares polynomial (of dimension U) epigraph of ℓ₁ norm (of dimension R) cone, parametrized by vector of matrices Ps derived from interpolant basis and polynomial domain constraints.

WSOSInterpEpiNormOne{T}(R::Int, U::Int, Ps::Vector{Matrix{T}}, use_dual::Bool = false)
source

Helpers for PosSemidefTriSparse

EpiPerSepSpectral helpers

Hypatia.Cones.MatrixCSqrType
struct MatrixCSqr{T<:Real, R<:Union{Complex{T<:Real}, T<:Real}} <: Hypatia.Cones.ConeOfSquares{T<:Real}

Real symmetric or complex Hermitian positive semidefinite cone of squares.

source
Hypatia.Cones.vector_dimFunction
vector_dim(
    _::Type{var"#s68"} where var"#s68"<:Hypatia.Cones.VectorCSqr,
    d::Int64
) -> Int64

The rank of the vector cone of squares, equal to the vector length.

source
vector_dim(
    _::Type{var"#s67"} where var"#s67"<:(Hypatia.Cones.MatrixCSqr{var"#s66", R} where var"#s66"<:Real),
    d::Int64
) -> Int64

The rank of the matrix cone of squares, equal to the side dimension of the matrix.

source
Hypatia.Cones.NegSqrtSSFType
struct NegSqrtSSF <: Hypatia.Cones.SepSpectralFun

The negative square root function $x \to -x^{1/2}$. Note this is a special case of the negative power: NegPower01SSF(0.5).

source
Hypatia.Cones.NegPower01SSFType
struct NegPower01SSF <: Hypatia.Cones.SepSpectralFun

The negative power function $x \to -x^p$ parametrized by $p \in (0, 1)$

source
Hypatia.Cones.Power12SSFType
struct Power12SSF <: Hypatia.Cones.SepSpectralFun

The power function $x \to x^p$ parametrized by $p \in (1, 2]$. Note for $p = 2$, it is more efficient to use EpiPerSquare.

source