Constraints

Equality constraints between polynomials

Equality between polynomials in PolyJuMP uses the same syntax as equality between affine or quadratic expression in JuMP. For instance, creating two quadratic n-variate polynomials p and q that must sum up to one can be done as follows:

julia> n = 3
3

julia> using DynamicPolynomials

julia> @polyvar x[1:n]
(Variable{DynamicPolynomials.Commutative{DynamicPolynomials.CreationOrder}, Graded{LexOrder}}[x₁, x₂, x₃],)

julia> X = monomials(x, 0:2)
10-element MonomialVector{DynamicPolynomials.Commutative{DynamicPolynomials.CreationOrder}, Graded{LexOrder}}:
 1
 x₃
 x₂
 x₁
 x₃²
 x₂x₃
 x₂²
 x₁x₃
 x₁x₂
 x₁²

julia> using SumOfSquares

julia> model = Model()
A JuMP Model
Feasibility problem with:
Variables: 0
Model mode: AUTOMATIC
CachingOptimizer state: NO_OPTIMIZER
Solver name: No optimizer attached.

julia> @variable(model, p, Poly(X))
(_[1]) + (_[2])x₃ + (_[3])x₂ + (_[4])x₁ + (_[5])x₃² + (_[6])x₂x₃ + (_[7])x₂² + (_[8])x₁x₃ + (_[9])x₁x₂ + (_[10])x₁²

julia> @variable(model, q, Poly(X))
(_[11]) + (_[12])x₃ + (_[13])x₂ + (_[14])x₁ + (_[15])x₃² + (_[16])x₂x₃ + (_[17])x₂² + (_[18])x₁x₃ + (_[19])x₁x₂ + (_[20])x₁²

julia> @constraint(model, p + q == 1)
(_[1] + _[11] - 1) + (_[2] + _[12])x₃ + (_[3] + _[13])x₂ + (_[4] + _[14])x₁ + (_[5] + _[15])x₃² + (_[6] + _[16])x₂x₃ + (_[7] + _[17])x₂² + (_[8] + _[18])x₁x₃ + (_[9] + _[19])x₁x₂ + (_[10] + _[20])x₁² ∈ PolyJuMP.ZeroPoly()

Vectorized constraints can also be used as well as vector of constraints, named constraints, ... For instance, if P and Q are two $n \times n$ matrices of polynomials, the following constraints the sum of rows and columns to match:

@constraint(model, con[i=1:n], sum(P[i, :]) == sum(Q[:, i]))

and con[i] contains the reference to the constraint between the ith row of P and the ith column of Q.

Inequality constraints between polynomials

Polynomials can be constrained to be sum-of-squares with the in syntax. For instance, to constrain a polynomial p to be sum-of-squares, do

julia> @constraint(model, p in SOSCone())
(_[1]) + (_[2])x₃ + (_[3])x₂ + (_[4])x₁ + (_[5])x₃² + (_[6])x₂x₃ + (_[7])x₂² + (_[8])x₁x₃ + (_[9])x₁x₂ + (_[10])x₁² is SOS

Automatically interpreting polynomial nonnegativity as a sum-of-squares constraint

As detailed in When is nonnegativity equivalent to sum of squares ?, the nonnegativity of a polynomial is not equivalent to the existence of a sum-of-squares decomposition. However, if explicitely specified, nonnegativity constraints can be automatically interpreted as sum-of-squares constraints. The simplest way to do that is to create the model with

model = SOSModel(...)

instead of

model = Model(...)

An alternative equivalent way is to call setpolymodule! after creating the model:

julia> setpolymodule!(model, SumOfSquares)

This second approach may be useful if the SumOfSquares JuMP extension need to be used with another JuMP extension that also has a special model constructor. A third alternative is the following:

julia> PolyJuMP.setdefault!(model, PolyJuMP.NonNegPoly, SOSCone)
SOSCone (alias for NonnegPolyInnerCone{MathOptInterface.PositiveSemidefiniteConeTriangle})

julia> PolyJuMP.setdefault!(model, PolyJuMP.PosDefPolyMatrix, SOSMatrixCone)
SOSMatrixCone (alias for PSDMatrixInnerCone{MathOptInterface.PositiveSemidefiniteConeTriangle})

This approach adds the flexibility to choose the default cone for

  • constraints of the form @constraint(mode, ..., some_polynomial ≥ other_polynomial, ...) which is the cone given as default to PolyJuMP.NonNegPoly; and
  • constraints of the form @constraint(mode, ..., some_matrix_of_polynomial in PSDCone(), ...) or @constraint(mode, ..., some_matrix_of_polynomial >= other_matrix_of_polynomial, PSDCone(), ...) which is the cone given as default to PolyJuMP.NonNegPolyMatrix.

For instance, to use the diagonally-dominant-sum-of-squares cone (see [Definition 2, AM17]) for the first type of contraints, do

julia> PolyJuMP.setdefault!(model, PolyJuMP.NonNegPoly, DSOSCone)
DSOSCone (alias for NonnegPolyInnerCone{SumOfSquares.DiagonallyDominantConeTriangle})

Changing the polynomial basis

As introduced in Choosing a polynomial basis, there may be numerical advantages to use another basis than the standard monomial basis when creating polynomial variables. Similarly, other polynomial bases can be used for polynomial constraints. However, for constraints, the polynomial space is determined by the polynomial constrained to be nonnegative. For instance, consider the constraint:

julia> using DynamicPolynomials

julia> @polyvar x y
(x, y)

julia> using SumOfSquares

julia> model = SOSModel()
A JuMP Model
Feasibility problem with:
Variables: 0
Model mode: AUTOMATIC
CachingOptimizer state: NO_OPTIMIZER
Solver name: No optimizer attached.

julia> @variable(model, α)
α

julia> @variable(model, β)
β

julia> @constraint(model, α * x^2 + β * y^2 ≥ (α - β) * x * y)
(β)y² + (-α + β)xy + (α)x² is SOS

where α and β are JuMP decision variables and x and y are polynomial variables. Since the polynomial is a quadratic form, the sum-of-squares certificate is also a quadratic form (see [Section~3.3.4, BPT12]). Hence the default polynomial basis used for the [Nonnegative polynomial variables] certificate is MonomialBasis([x, y]), that is, we search for a positive semidefinite matrix Q such that

\[\alpha x^2 + \beta y^2 - (\alpha - \beta) x y = X^\top Q X\]

where $X = (x, y)$.

As the polynomial space is determined by the polynomial being constrained, only the basis type needs to be given. For instance, to use the scaled monomial basis in the example above, use

julia> @constraint(model, α * x^2 + β * y^2 ≥ (α - β) * x * y, basis = ScaledMonomialBasis)
(β)y² + (-α + β)xy + (α)x² is SOS

Polynomial nonnegativity on a subset of the space

By default, the constraint

julia> @constraint(model, x^3 - x^2 + 2x*y -y^2 + y^3 >= α)
(-α) + (-1)y² + (2)xy + (-1)x² + (1)y³ + (1)x³ is SOS

constrains the polynomial to be nonnegative for every real numbers x and y. However, the set of points (x, y) for which the polynomial is constrained to be nonnegative can be specified by the domain keyword:

julia> S = @set x >= 0 && y >= 0 && x + y >= 1;

julia> @constraint(model, x^3 - x^2 + 2x*y -y^2 + y^3 >= α, domain = S)
(-α) + (-1)y² + (2)xy + (-1)x² + (1)y³ + (1)x³ is SOS

See this notebook for a detailed example.

Dual of polynomial constraints

The dual of a polynomial constraint cref is a moment serie μ as defined in MultivariateMoments. The dual can be obtained with the dual function as with classical dual values in JuMP.

μ = dual(cref)

By dual of a Sum-of-Squares constraint, we may mean different things and the meaning chosen for dual function was chosen for consistency with the definition of the JuMP dual function to ensure that generic code will work as expected with Sum-of-Squares constraints. In a Sum-of-Squares constraint, a polynomial $p$ is constraint to be SOS in some domain defined by polynomial q_i. So p(x) is constrained to be equal to s(x) = s_0(x) + s_1(x) * q_1(x) + s_2(x) * q_2(x) + ... where the s_i(x) polynomials are Sum-of-Squares. The dual of the equality constraint between p(x) and s(x) is given by SumOfSquares.MultivariateMoments.moments.

μ = moments(cref)

Note that the dual and moments may give different results. For instance, the output of dual only contains the moments corresponding to monomials of p while the output of moments may give the moments of other monomials if s(x) has more monomials than p(x). Besides, if the domain contains polynomial, equalities, only the remainder of p(x) - s(x) modulo the ideal is constrained to be zero, see Corollary 2 of [CLO13]. In that case, the output moments is the dual of the constraint on the remainder so some monomials may have different moments with dual or moments.

The dual of the Sum-of-Squares constraint on s_0(x), commonly referred to as the the matrix of moments can be obtained using moment_matrix:

ν = moment_matrix(cref)

The atomic_measure function of MultivariateMoments can be used to check if there exists an atomic measure (i.e. a measure that is a sum of Dirac measures) that has the moments given in the the moment matrix ν. This can be used for instance in polynomial optimization (see this notebook) or stability analysis (see this notebook).

References

[BPT12] Blekherman, G.; Parrilo, P. A. & Thomas, R. R. Semidefinite Optimization and Convex Algebraic Geometry. Society for Industrial and Applied Mathematics, 2012.

[CLO13] Cox, D., Little, J., & OShea, D. Ideals, varieties, and algorithms: an introduction to computational algebraic geometry and commutative algebra. Springer Science & Business Media, 2013.

[AM17] Ahmadi, A. A. & Majumdar, A. DSOS and SDSOS Optimization: More Tractable Alternatives to Sum of Squares and Semidefinite Optimization. ArXiv e-prints, 2017.

API Reference

Default choice for the maxdegree keyword:

SumOfSquares.default_maxdegreeFunction
default_maxdegree(monos, domain)

Return the default maxdegree to use for certifying a polynomial with monomials monos to be Sum-of-Squares over the domain domain. By default, the maximum of the maxdegree of monos and of all multipliers of the domain are used so that at least constant multipliers can be used with a Putinar certificate.

source

Special case that is second-order cone representable:

Inner approximations of the PSD cone that do not require semidefinite programming:

SumOfSquares.DiagonallyDominantConeTriangleType
struct DiagonallyDominantConeTriangle <: MOI.AbstractSymmetricMatrixSetTriangle
    side_dimension::Int
end

See Definition 4 of [AM17] for a precise definition of the last two items.

[AM17] Ahmadi, A. A. & Majumdar, A. DSOS and SDSOS Optimization: More Tractable Alternatives to Sum of Squares and Semidefinite Optimization ArXiv e-prints, 2017.

source
SumOfSquares.ScaledDiagonallyDominantConeTriangleType
struct ScaledDiagonallyDominantConeTriangle <: MOI.AbstractSymmetricMatrixSetTriangle
    side_dimension::Int
end

See Definition 4 of [AM17] for a precise definition of the last two items.

[AM17] Ahmadi, A. A. & Majumdar, A. DSOS and SDSOS Optimization: More Tractable Alternatives to Sum of Squares and Semidefinite Optimization ArXiv e-prints, 2017.

source

Approximations of the cone of nonnegative polynomials:

SumOfSquares.NonnegPolyInnerConeType
struct NonnegPolyInnerCone{MCT <: MOI.AbstractVectorSet}
end

Inner approximation of the cone of nonnegative polynomials defined by the set of polynomials of the form

X^\top Q X

where X is a vector of monomials and Q is a symmetric matrix that belongs to the cone psd_inner.

source
SumOfSquares.SDSOSConeType
const SDSOSCone = NonnegPolyInnerCone{ScaledDiagonallyDominantConeTriangle}

Scaled-diagonally-dominant-sum-of-squares cone; see [Definition 2, AM17] and NonnegPolyInnerCone.

[AM17] Ahmadi, A. A. & Majumdar, A. DSOS and SDSOS Optimization: More Tractable Alternatives to Sum of Squares and Semidefinite Optimization ArXiv e-prints, 2017.

source
SumOfSquares.DSOSConeType
const DSOSCone = NonnegPolyInnerCone{DiagonallyDominantConeTriangle}

Diagonally-dominant-sum-of-squares cone; see [Definition 2, AM17] and NonnegPolyInnerCone.

[AM17] Ahmadi, A. A. & Majumdar, A. DSOS and SDSOS Optimization: More Tractable Alternatives to Sum of Squares and Semidefinite Optimization ArXiv e-prints, 2017.

source

Approximations of the cone of positive semidefinite polynomial matrices:

SumOfSquares.PSDMatrixInnerConeType
struct PSDMatrixInnerCone{MCT <: MOI.AbstractVectorSet} <: PolyJuMP.PolynomialSet
end

Inner approximation of the cone of polynomial matrices P(x) that are positive semidefinite for any x defined by the set of $n \times n$ polynomial matrices such that the polynomial $y^\top P(x) y$ belongs to NonnegPolyInnerCone{MCT} where y is a vector of $n$ auxiliary polynomial variables.

source
SumOfSquares.SOSMatrixConeType
const SOSMatrixCone = PSDMatrixInnerCone{MOI.PositiveSemidefiniteConeTriangle}

Sum-of-squares matrices cone; see [Section 3.3.2, BPT12] and PSDMatrixInnerCone.

[BPT12] Blekherman, G.; Parrilo, P. A. & Thomas, R. R. Semidefinite Optimization and Convex Algebraic Geometry. Society for Industrial and Applied Mathematics, 2012.

source

Approximations of the cone of convex polynomials:

SumOfSquares.ConvexPolyInnerConeType
struct ConvexPolyInnerCone{MCT} <: PolyJuMP.PolynomialSet end

Inner approximation of the set of convex polynomials defined by the set of polynomials such that their hessian belongs to to the set PSDMatrixInnerCone{MCT}()

source
SumOfSquares.SOSConvexConeType
const SOSConvexCone = ConvexPolyInnerCone{MOI.PositiveSemidefiniteConeTriangle}

Sum-of-squares convex polynomials cone; see [Section 3.3.3, BPT12] and ConvexPolyInnerCone.

[BPT12] Blekherman, G.; Parrilo, P. A. & Thomas, R. R. Semidefinite Optimization and Convex Algebraic Geometry. Society for Industrial and Applied Mathematics, 2012.

source

Approximations of the cone of copositive matrices:

SumOfSquares.CopositiveInnerType
struct CopositiveInner{S} <: PolyJuMP.PolynomialSet
    # Inner approximation of the PSD cone, i.e. typically either
    # `SOSCone`, `DSOSCone` or `SDSOSCone`,
    psd_inner::S
end

A symmetric matrix $Q$ is copositive if $x^{\top} Q x \ge 0$ for all vector $x$ in the nonnegative orthant. Checking copositivity is a co-NP-complete problem [MK87] and this cone is only the inner approximation of the cone of copositive symmetric matrices given by Minknowski sum of psd_inner and the cone of symmetric matrices with nonnegative entries (the diagonal entries can be chosen to be zero) [Lemma 3.164, BPT12].

The matrix with nonnegative entries can be interpreted as lagrangian multipliers. For instance,

@polyvar x y
@constraint(model, x^2 - 2x*y + y^2 in CopositiveInner(SOSCone()))

is equivalent to

# Matrix that we require to be copositive
Q = [ 1 -1
     -1  1]
λ = @variable(model, lower_bound=0)
# Symmetric matrix of nonnegative entries
Λ = [0 λ
     λ 0]
using LinearAlgebra # For `Symmetric`
@constraint(model, Symmetric(Q - Λ) in PSDCone())

which is equivalent to

@polyvar x y
λ = @variable(model, lower_bound=0)
@constraint(model, x^2 - 2x*y + y^2 - 2*λ * x*y in SOSCone())

which is the same as, using the domain keyword,

@polyvar x y
@constraint(model, x^2 - 2x*y + y^2 in SOSCone(), domain = @set x*y ≥ 0)

As an important difference with its equivalent forms, the GramMatrixAttribute for the copositive constraint is given by matrix Q while for the equivalent form using the domainkeyword, the value of the attribute would correspond to the the gram matrix in thepsd_innercone, i.e. which should be equal toQ - Λ`.

[BPT12] Blekherman, G.; Parrilo, P. A. & Thomas, R. R. Semidefinite Optimization and Convex Algebraic Geometry. Society for Industrial and Applied Mathematics, 2012.

[MK87] K. G. Murty and S. N. Kabadi. Some NP-complete problems in quadratic and nonlinear programming. Mathematical programming, 39:117–129, 1987.

source

Types of sparsity

SumOfSquares.Certificate.Sparsity.VariableType
struct Sparsity.Variable <: Sparsity.Pattern end

Variable or correlative sparsity as developed in [WSMM06].

[WSMM06] Waki, Hayato, Sunyoung Kim, Masakazu Kojima, and Masakazu Muramatsu. "Sums of squares and semidefinite program relaxations for polynomial optimization problems with structured sparsity." SIAM Journal on Optimization 17, no. 1 (2006): 218-242.

source
SumOfSquares.Certificate.Sparsity.MonomialType
struct Sparsity.Monomial{C<:CEG.AbstractCompletion} <: Sparsity.Pattern
    completion::C
    k::Int
    use_all_monomials::Bool
end

Monomial or term sparsity as developed in [WML20a, WML20b]. The completion field should be ClusterCompletion() [default] for the block-closure or cluster completion [WML20a], and ChordalCompletion() for chordal completion [WML20b]. The integer k [default=0] corresponds to Σ_k defined in [(3.2), WML20a] and k = 0 corresponds to Σ_* defined in [(3.3), WML20a]. If use_all_monomials is false then some monomials of the basis might be dropped from the basis if not needed.

[WML20a] Wang, Jie, Victor Magron, and Jean-Bernard Lasserre. TSSOS: A Moment-SOS hierarchy that exploits term sparsity. arXiv preprint arXiv:1912.08899 (2020).

[WML20b] Wang, Jie, Victor Magron, and Jean-Bernard Lasserre. Chordal-TSSOS: a moment-SOS hierarchy that exploits term sparsity with chordal extension. arXiv preprint arXiv:2003.03210 (2020).

source
SumOfSquares.Certificate.Sparsity.SignSymmetryType
struct SignSymmetry <: Sparsity.Pattern end

Sign symmetry as developed in [Section III.C, L09]. Let n be the number of variables. A sign-symmetry is a binary vector r of {0, 1}^n such that dot(r, e) is even for all exponent e.

Let o(e) be the binary vector of {0, 1}^n for which the ith bit is i iff the ith entry of e is odd. Let O be the set of o(e) for exponents of e. A binary vector r is a sign-symmetry if an only if it is orthogonal to all the elements of O.

Since we are only interested in the orthogonal subspace, say R, of O, even if O is not a linear subspace (i.e., it is not invariant under xor), we compute its span. We start by creating row echelon form of the span of O using XORSpace. We then compute a basis for R. The set R of sign symmetries will be the span of this basis.

Let a, b be exponents of monomials of the gram basis. For any r in R,

⟨r, o(a + b)⟩ = ⟨r, xor(o(a), o(b)⟩ = xor(⟨r, o(a)⟩, ⟨r, o(b)⟩)

For o(a, b) to be sign symmetric, this scalar product should be zero for all sign symmetry r. This is equivalent to saying that ⟨r, o(a)⟩ and ⟨r, o(b)⟩ are equal for all r in R. In other words, the projection of o(a) and o(b) to R have the same coordinates in the basis. If we order the monomials by grouping them by equal coordinates of projection, we see that the product that are sign symmetric form a block diagonal structure. This means that we can group the monomials by these coordinates.

[L09] Löfberg, Johan. Pre-and post-processing sum-of-squares programs in practice. IEEE transactions on automatic control 54, no. 5 (2009): 1007-1011.

source
SumOfSquares.Certificate.Sparsity.XORSpaceType
mutable struct XORSpace{T<:Integer}
    dimension::Int
    basis::Vector{T}
end

Basis for a linear subspace of the Hamming space (i.e. set of binary string {0, 1}^n of length n) represented in the bits of an integer of type T. This is used to implement Certificate.Sparsity.SignSymmetry.

Consider the scalar product ⟨a, b⟩ which returns the the xor of the bits of a & b. It is a scalar product since ⟨a, b⟩ = ⟨b, a⟩ and ⟨a, xor(b, c)⟩ = xor(⟨a, b⟩, ⟨a, c⟩).

We have two options here to compute the orthogonal space.

The first one is to build an orthogonal basis with some kind of Gram-Schmidt process and then to obtain the orthogonal space by removing the projection from each vector of the basis.

The second option, which is the one we use here is to compute the row echelon form and then read out the orthogonal subspace directly from it. For instance, if the row echelon form is

1 a 0 c e
  b 1 d f

then the orthogonal basis is

a 1 b 0 0
c 0 d 1 0
e 0 f 0 1

The basis vector has dimension nonzero elements. Any element added with push! can be obtained as the xor of some of the elements of basis. Moreover, the basis is kept in row echelon form in the sense that the first, second, ..., i - 1th bits of basis[i] are zero and basis[i] is zero if its ith bit is not one.

source

Ideal certificates:

SumOfSquares.Certificate.MaxDegreeType
struct MaxDegree{CT <: SumOfSquares.SOSLikeCone, BT <: MB.AbstractPolynomialBasis} <: SimpleIdealCertificate{CT, BT}
    cone::CT
    basis::Type{BT}
    maxdegree::Int
end

The MaxDegree certificate ensures the nonnegativity of p(x) for all x such that h_i(x) = 0 by exhibiting a Sum-of-Squares polynomials σ(x) such that p(x) - σ(x) is guaranteed to be zero for all x such that h_i(x) = 0. The polynomial σ(x) is search over cone with a basis of type basis such that the degree of σ(x) does not exceed maxdegree.

source
SumOfSquares.Certificate.FixedBasisType
struct FixedBasis{CT <: SumOfSquares.SOSLikeCone, BT <: MB.AbstractPolynomialBasis} <: SimpleIdealCertificate{CT, BT}
    cone::CT
    basis::BT
end

The FixedBasis certificate ensures the nonnegativity of p(x) for all x such that h_i(x) = 0 by exhibiting a Sum-of-Squares polynomials σ(x) such that p(x) - σ(x) is guaranteed to be zero for all x such that h_i(x) = 0. The polynomial σ(x) is search over cone with basis basis.

source
SumOfSquares.Certificate.NewtonType
struct Newton{CT <: SumOfSquares.SOSLikeCone, BT <: MB.AbstractPolynomialBasis, NPT <: Tuple} <: SimpleIdealCertificate{CT, BT}
    cone::CT
    basis::Type{BT}
    variable_groups::NPT
end

The Newton certificate ensures the nonnegativity of p(x) for all x such that h_i(x) = 0 by exhibiting a Sum-of-Squares polynomials σ(x) such that p(x) - σ(x) is guaranteed to be zero for all x such that h_i(x) = 0. The polynomial σ(x) is search over cone with a basis of type basis chosen using the multipartite Newton polytope with parts variable_groups. If variable_groups = tuple() then it falls back to the classical Newton polytope with all variables in the same part.

source
SumOfSquares.Certificate.RemainderType
struct Remainder{GCT<:AbstractIdealCertificate} <: AbstractIdealCertificate
    gram_certificate::GCT
end

The Remainder certificate ensures the nonnegativity of p(x) for all x such that h_i(x) = 0 by guaranteeing the remainder of p(x) modulo the ideal generated by ⟨h_i⟩ to be nonnegative for all x such that h_i(x) = 0 using the certificate gram_certificate. For instance, if gram_certificate is SumOfSquares.Certificate.Newton, then the certificate Remainder(gram_certificate) will take the remainder before computing the Newton polytope hence might generate a much smaller Newton polytope hence a smaller basis and smaller semidefinite program. However, this then corresponds to a lower degree of the hierarchy which might be insufficient to find a certificate.

source
SumOfSquares.Certificate.Sparsity.IdealType
struct Sparsity.Ideal{S <: Sparsity.Pattern, C <: AbstractIdealCertificate} <: SumOfSquares.Certificate.AbstractIdealCertificate
    sparsity::S
    certificate::C
end

Same certificate as certificate except that the Sum-of-Squares polynomial σ is modelled as a sum of Sum-of-Squares polynomials with smaller bases using the sparsity reduction sparsity.

source
SumOfSquares.Certificate.Symmetry.IdealType
struct Symmetry.Ideal{C,GT,AT<:SymbolicWedderburn.Action} <: SumOfSquares.Certificate.AbstractIdealCertificate
    pattern::Symmetry.Pattern{GT,AT}
    certificate::C
end

Same certificate as certificate except that the Sum-of-Squares polynomial σ is modelled as a sum of Sum-of-Squares polynomials with smaller bases using the Symbolic Wedderburn decomposition of the symmetry pattern specified by pattern.

source

Preorder certificates:

SumOfSquares.Certificate.PutinarType
struct Putinar{
    MC<:AbstractIdealCertificate,
    IC<:AbstractIdealCertificate,
} <: AbstractPreorderCertificate
    multipliers_certificate::MC
    ideal_certificate::IC
    maxdegree::Int
end

The Putinar certificate ensures the nonnegativity of p(x) for all x such that g_i(x) >= 0 and h_i(x) = 0 by exhibiting Sum-of-Squares polynomials σ_i(x) such that p(x) - ∑ σ_i(x) g_i(x) is guaranteed to be nonnegativity for all x such that h_i(x) = 0. The polynomials σ_i(x) are search over cone with a basis of type basis such that the degree of σ_i(x) g_i(x) does not exceed maxdegree.

source
SumOfSquares.Certificate.Sparsity.PreorderType
struct Sparsity.Preorder{S <: Sparsity.Pattern, C <: SumOfSquares.Certificate.AbstractPreorderCertificate} <: SumOfSquares.Certificate.AbstractPreorderCertificate
    sparsity::S
    certificate::C
end

Same certificate as C except that the Sum-of-Squares polynomials σ_i are modelled as a sum of Sum-of-Squares polynomials with smaller basis using the sparsity reduction sparsity.

source

Attributes

PolyJuMP.MomentsAttributeType
MomentsAttribute(N)
MomentsAttribute()

A constraint attribute for the vector of moments corresponding to the constraint that a polynomial is zero in the full space in result N. If the polynomial is constrained to be zero in an algebraic set, it is the moments for the constraint once it is rewritten into an constraint on the full space. If N is omitted, it is 1 by default.

Examples

Consider the following program:

@variable(model, α)
@variable(model, β ≤ 1)

using DynamicPolynomials
@polyvar x y
cref = @constraint(model, α * x - β * y == 0, domain = @set x == y)

The constraint is equivalent to

@constraint(model, (α - β) * y == 0)

for which the dual is the 1-element vector with the moment of y of value -1. This is the result of moments(cref). However, the dual of cref obtained by dual(cref) is the 2-elements vector with both the moments of x and y of value -1.

source
SumOfSquares.GramMatrixType
struct GramMatrix{T, B, U, MT <: AbstractMatrix{T}} <: AbstractGramMatrix{T, B, U}
    Q::MT
    basis::B
end

Gram matrix $x^\top Q x$ where Q is a symmetric matrix indexed by the vector of polynomials of the basis basis.

source
SumOfSquares.GramMatrixAttributeType
GramMatrixAttribute(N)
GramMatrixAttribute()

A constraint attribute for the GramMatrix of a constraint, that is, the positive semidefinite matrix Q indexed by the monomials in the vector X such that $X^\top Q X$ is the sum-of-squares certificate of the constraint.

source
SumOfSquares.gram_operateFunction
gram_operate(::typeof(+), p::GramMatrix, q::GramMatrix)

Computes the Gram matrix equal to the sum between p and q. On the opposite, p + q gives a polynomial equal to p + q. The polynomial p + q can also be obtained by polynomial(gram_operate(+, p, q)).

source
gram_operate(/, p::GramMatrix, α)

Computes the Gram matrix equal to p / α. On the opposite, p / α gives a polynomial equal to p / α. The polynomial p / α can also be obtained by polynomial(gram_operate(/, p, α)).

source
SumOfSquares.LagrangianMultipliersType
LagrangianMultipliers(N)
LagrangianMultipliers()

A constraint attribute fot the LagrangianMultipliers associated to the inequalities of the domain of a constraint. There is one multiplier per inequality and no multiplier for equalities as the equalities are handled by reducing the polynomials over the ideal they generate instead of explicitely creating multipliers.

source

Bridges are automatically added using the following utilities:

PolyJuMP.bridgeableFunction
bridgeable(c::JuMP.AbstractConstraint, S::Type{<:MOI.AbstractSet})

Wrap the constraint c in JuMP.BridgeableConstraints that may be needed to bridge variables constrained in S on creation.

bridgeable(c::JuMP.AbstractConstraint, F::Type{<:MOI.AbstractFunction},
           S::Type{<:MOI.AbstractSet})

Wrap the constraint c in JuMP.BridgeableConstraints that may be needed to bridge F-in-S constraints.

source
PolyJuMP.bridgesFunction
bridges(F::Type{<:MOI.AbstractFunction}, S::Type{<:MOI.AbstractSet})

Return a list of bridges that may be needed to bridge F-in-S constraints but not the bridges that may be needed by constraints added by the bridges.

source
bridges(S::Type{<:MOI.AbstractSet})

Return a list of bridges that may be needed to bridge variables constrained in S on creation but not the bridges that may be needed by constraints added by the bridges.

source

Chordal extension:

SumOfSquares.Certificate.Sparsity.ChordalExtensionGraph.completionFunction
completion(G::Graph, comp::ChordalCompletion)

Return a chordal extension of G and the corresponding maximal cliques.

The algoritm is Algorithm 3 in [BA10] with the GreedyFillIn heuristic of Table I.

[BA10] Bodlaender, Hans L., and Arie MCA Koster. Treewidth computations I. Upper bounds. Information and Computation 208, no. 3 (2010): 259-275. Utrecht University, Utrecht, The Netherlands www.cs.uu.nl

source
completion(G::Graph, comp::ChordalCompletion)

Return a cluster completion of G and the corresponding maximal cliques.

source

Bridges for polynomial optimization

PolyJuMP.ScalarPolynomialFunctionType
struct ScalarPolynomialFunction{T,P<:MP.AbstractPolynomial{T}} <: MOI.AbstractScalarFunction
    polynomial::P
    variables::Vector{MOI.VariableIndex}
end

Defines the polynomial function of the variables variables where the variable variables(p)[i] corresponds to variables[i].

source
PolyJuMP.Bridges.Objective.ToPolynomialBridgeType
ToPolynomialBridge{T}

ToPolynomialBridge implements the following reformulations:

  • $\min \{f(x)\}$ into $\min\{p(x)\}$
  • $\max \{f(x)\}$ into $\max\{p(x)\}$

where $f(x)$ is a scalar function and $p(x)$ is a polynomial.

Source node

ToPolynomialBridge supports:

  • MOI.ObjectiveFunction{F} where F is a MOI.AbstractScalarFunction for which convert(::Type{PolyJuMP.ScalarPolynomialFunction}, ::Type{F}). That is for instance the case for MOI.VariableIndex, MOI.ScalarAffineFunction and MOI.ScalarQuadraticFunction.

Target nodes

ToPolynomialBridge creates:

  • One objective node: MOI.ObjectiveFunction{PolyJuMP.ScalarPolynomialFunction{T}}
source
PolyJuMP.Bridges.Constraint.ToPolynomialBridgeType
ToPolynomialBridge{T,S} <: Bridges.Constraint.AbstractBridge

ToPolynomialBridge implements the following reformulations:

  • $f(x) \in S$ into $p(x) \in S$ where $f(x)$ is a scalar function and $p(x)$ is a polynomial.

Source node

ToPolynomialBridge supports:

  • F in S where F is a MOI.AbstractScalarFunction for which

convert(::Type{PolyJuMP.ScalarPolynomialFunction}, ::Type{F}). That is for instance the case for MOI.VariableIndex, MOI.ScalarAffineFunction and MOI.ScalarQuadraticFunction.

Target nodes

ToPolynomialBridge creates:

source

SAGE extension

To use the SAGE cone in place of the Sum-of-Squares cone for an inequality constraints between polynomials, use the following:

import PolyJuMP
PolyJuMP.setpolymodule!(model, PolyJuMP.SAGE)
PolyJuMP.SAGE.PolynomialsType
struct Polynomials{M<:Union{Nothing,Int,MP.AbstractMonomial}} <: PolyJuMP.PolynomialSet

Sums of AM/GM Exponential for polynomials.

source
PolyJuMP.SAGE.SignomialsType
struct Signomials{M<:Union{Nothing,Int,MP.AbstractMonomial}} <: PolyJuMP.PolynomialSet

Sums of AM/GM Exponential for signomials.

source
PolyJuMP.SAGE.SignomialsBridgeType
SignomialsBridge{T,S,P,F} <: MOI.Bridges.Constraint.AbstractBridge

We use the Signomials Representative SR equation of [MCW21].

[MCW20] Riley Murray, Venkat Chandrasekaran, Adam Wierman "Newton Polytopes and Relative Entropy Optimization" https://arxiv.org/abs/1810.01614 [MCW21] Murray, Riley, Venkat Chandrasekaran, and Adam Wierman. "Signomials and polynomial optimization via relative entropy and partial dualization." Mathematical Programming Computation 13 (2021): 257-295. https://arxiv.org/pdf/1907.00814.pdf

source
PolyJuMP.SAGE.AGEBridgeType
AGEBridge{T,F,G,H} <: MOI.Bridges.Constraint.AbstractBridge

The nonnegativity of x ≥ 0 in

∑ ci x^αi ≥ -c0 x^α0

can be reformulated as

∑ ci exp(αi'y) ≥ -β exp(α0'y)

In any case, it is shown to be equivalent to

∃ ν ≥ 0 s.t. D(ν, exp(1)*c) ≤ β, ∑ νi αi = α0 ∑ νi [CP16, (3)]

where N(ν, λ) = ∑ νj log(νj/λj) is the relative entropy function. The constant exp(1) can also be moved out of D into

∃ ν ≥ 0 s.t. D(ν, c) - ∑ νi ≤ β, ∑ νi αi = α0 ∑ νi [MCW21, (2)]

The relative entropy cone in MOI is (u, v, w) such that D(w, v) ≥ u. Therefore, we can either encode (β, exp(1)*c, ν) [CP16, (3)] or (β + ∑ νi, c, ν) [MCW21, (2)]. In this bridge, we use the second formulation.

Note

A direct application of the Arithmetic-Geometric mean inequality gives

∃ ν ≥ 0 s.t. D(ν, exp(1)*c) ≤ -log(-β), ∑ νi αi = α0, ∑ νi = 1 [CP16, (4)]

which is not jointly convex in (ν, c, β). The key to get the convex formulation [CP16, (3)] or [MCW21, (2)] instead is to use the convex conjugacy between the exponential and the negative entropy functions [CP16, (iv)].

[CP16] Chandrasekaran, Venkat, and Parikshit Shah. "Relative entropy relaxations for signomial optimization." SIAM Journal on Optimization 26.2 (2016): 1147-1173. [MCW21] Murray, Riley, Venkat Chandrasekaran, and Adam Wierman. "Signomials and polynomial optimization via relative entropy and partial dualization." Mathematical Programming Computation 13 (2021): 257-295. https://arxiv.org/pdf/1907.00814.pdf

source

Internal functions

SumOfSquares.Certificate.Symmetry._reorder!Function
_reorder!(F::LinearAlgebra.Schur{T}) where {T}

Given a Schur decomposition of a, reorder it so that its eigenvalues are in in increasing order.

Note that if T<:Real, F.Schur is quasi upper triangular. By (quasi), we mean that there may be nonzero entries in S[i+1,i] representing complex conjugates. In that case, the complex conjugate are permuted together. If T<:Complex, then S is triangular.

source
SumOfSquares.Certificate.Symmetry._rotate_complexFunction
_rotate_complex(A::AbstractMatrix{T}, B::AbstractMatrix{T}; tol = Base.rtoldefault(real(T))) where {T}

Given (quasi) upper triangular matrix A and B that have the eigenvalues in the same order except the complex pairs which may need to be (signed) permuted, returns an othogonal matrix P such that P' * A * P and B have matching low triangular part. The upper triangular part will be dealt with by _sign_diag.

By (quasi), we mean that if S is a Matrix{<:Real}, then there may be nonzero entries in S[i+1,i] representing complex conjugates. If S is a Matrix{<:Complex}, then S is upper triangular so there is nothing to do.

source