Backend

PolynomialOptimization mostly uses external solvers that have Julia bindings or are implemented in Julia, but it also provides own solver implementations particularly for the purpose of polynomial optimization. This page is only relevant if you intend to implement an interface between PolynomialOptimization and a new solver or if you want to provide missing functionality for an existing solver. It is of no relevance if you only want to use the existing solvers.

Warning

The package does not introduce any hard dependencies on the external solvers. Therefore, you may or may not decide to install them on your system. Instead, the solvers are addressed as weak dependencies. This means that you have to load the solver Julia package manually before you are able to use it for optimization.

PolynomialOptimization also provides an interface that can be implemented if another solver should be supported. This consists of just a few methods for the various functions.

MathOptInterface

Why do we re-invent the wheel and don't just use MathOptInterface? This would immediately bring compatibility with a huge number of solvers that are available under Julia.

This would certainly be possible. It might even be a not-too-bad idea to provide an automatic fallback to MathOptInterface in the future for solvers that are not supported natively.

However, MathOptInterface is enormous in its feature range. Considering that, it is amazingly fast. But in PolynomialOptimization, only a very small subset of all the features is required, and a lot of additional assumptions on how the problem is constructed can be made. This allows to use the solver's API in a most efficient way, which typically is not the way in which implementations of MathOptInterface would address the solver.

Additionally, in some situations, a lot of very similar sub-problems need to be solved; in these cases, PolynomialOptimization's own interface allows to keep the optimizer task alive and just do the tiny modification instead of setting things up from the start again - which for MathOptInterface only works if a solver is implemented in a particular way.

This focus on efficiency is very important as (relevant) polynomial optimization problems are huge. It is just not possible to waste time and memory in bookkeeping that is not really needed.

poly_optimize

The optimization of polynomial problems requires a solver that understands linear and semidefinite constraints. All functions in this section are defined (and exported) in the submodule PolynomialOptimization.Solver.

Solver interface

In general, a solver implementation can do whatever it wants; it just needs to implement the poly_optimize method with the appropriate Val-wrapped solver method as its first parameter. However, it is very helpful to just do some basic setup such as creating the solver object in this function and delegate all the work of setting up the actual problem to moment_setup! or sos_setup!. In order to do so, a solver implementation should create a new type that contains all the relevant data during setup of the problem. Usually, a solver falls in one of three categories:

  • Problem data has to be supplied in matrix/vector form; in this case, the new type should be a descendant of AbstractSparseMatrixSolver. Usually, open-source solvers fall in this category.
  • Problem data is constructed incrementally via various calls to API functions of the solver, which does not provide access to its internals. In this case, the new type should be a descendant of AbstractAPISolver. Usually, commercial solvers fall in this category.

However, it is not required that the type is in fact a subtype of either of those; the most general possible supertype is AbstractSolver, which does not make any assumptions or provide any but the most skeleton fallback implementations and a default for mindex.

PolynomialOptimization.Solver.AbstractSolverType
AbstractSolver{T,V<:Real}

Abstract supertype for any solver state. T is the type that is returned by calling mindex with such a state; V is the type of the coefficients used in the solver. Using the default implementation for mindex, T will be UInt. Most likely, V will be Float64.

source
PolynomialOptimization.Solver.mindexFunction
mindex(::AbstractSolver{T}, monomials::IntMonomialOrConj...)::T

Calculates the index that the product of all monomials will have in the SDP represented by state. The default implementation calculates the one-based monomial index according to a dense deglex order and returns an UInt. The returned index is arbitrary as long as it is unique for the total monomial.

source

Every implementation of poly_optimize should return a tuple that contains some internal state of the solver as well as the optimal value and the status of the solver. A method for issuccess should then translate this status into a simple boolean, where deciding on ambiguities (near success) is up to the solver implementation.

PolynomialOptimization.Solver.poly_optimizeMethod
poly_optimize(::Val{method}, relaxation::AbstractRelaxation,
    groupings::RelaxationGroupings; representation, verbose, kwargs...)

This is the central entry point that a solver has to implement. It has to carry out the optimization and must return a tuple of three values:

  1. An internal state that can be used later on to access the solver again to extract solutions or (if possible) reoptimization.
  2. The success status as returned by the solver and understood by the issuccess implementation.
  3. The minimum value of the objective.

See also moment_setup!, sos_setup!.

source
LinearAlgebra.issuccessMethod
issuccess(::Val{method}, status)

A solver must implement this method for all of its possible methods to indicate whether a status status signifies success.

source

Once a solver has been implemented, it should add its solver symbol to the vector solver_methods, which enables this solver to be chosen automatically. Apart from the exact specification :<solvername>Moment or :<solvername>SOS, a short form :<solvername> that chooses the recommended method should also be implemented. For this, the @solver_alias macro can be used. When details on the solution data a requested, the extract_moments, extract_sos, or extract_info function is called, where at least the former two have to be implemented for each solver:

PolynomialOptimization.Solver.extract_sosFunction
extract_sos(relaxation::AbstractRelaxation, state, ::Val{type},
    index::Union{<:Integer,<:AbstractUnitRange}, rawstate) where {type}

Extracts data that contains the raw solver information about the SOS data contained in the result. For moment optimizations, this corresponds to the dual data; for SOS optimizations, this is the primal data. rawstate is the return value of the preceding call to extract_sos_prepare (by default, nothing). Note that the SOS data may be queried in any order, partially or completely.

The parameters type and index indicates which constraint/variable the data corresponds to. type is a symbol, index is the range of indices within constraints of the same type, although both the type as well as the interpretation of index may change by providing custom definitions for addtocounter! or using the macros @counter_atomic and @counter_alias.

The return value of this function should be a scalar, vector, or matrix, depending on what data was requested. The following relations should hold:

typeresult type
:fix (moment only)vector
:free (SOS only)vector
:nonnegativevector
:quadraticvector
:rotated_quadraticvector
:psdvector[1] or matrix
:psd_complexvector[1][2] or matrix
:ddvector[1] or matrix
:dd_complexvector[1][2] or matrix
:lnormvector
:lnorm_complexvector[2]
:sddvector
:sdd_complexvector[2]
Info

It is guaranteed that the range that is queries using index always corresponds to data that was added contiguously, with no other cones interspersed.

source
PolynomialOptimization.Solver.extract_sos_prepareFunction
extract_sos_prepare(relaxation::AbstractRelaxation, state)

Prepares for one or multiple calls of extract_sos. The return value will be passed to the function as an argument. This is particularly relevant for solvers which don't allow the extraction of subsets of the data, but only the whole vector: Retrieve it here, then pass it as an output. The default implementation does nothing.

source

In order to relate aspects of the problem with data in the solver, the cones that are added are counted. This works automatically, keeping a separate counter for every type of cone and counting vector-valued cones (which are most) with their actual length. This behavior can be customized:

PolynomialOptimization.Solver.@counter_aliasMacro
@counter_alias(::Type{<:AbstractSolver}, counter, alias)

Defines the addtocounter! function in such a way that contraints of type counter instead only affect the counter alias. Do not use @counter_atomic together with this macro on the same counter.

counter may either be a Symbol, a tuple of Symbols, or the value Any. Note that Any has weakest precendence, irrespective of when the macro was called.

Warning

Regardless of whether counter or alias where made atomic before, after this macro, counter will not be so (although it shares the same counter as alias). This may or may not be desirable (most likely not), so always make atomic counters explicit using @counter_atomic, which allows the definition of aliases.

source
PolynomialOptimization.Solver.@counter_atomicMacro
@counter_atomic(::Type{<:AbstractSolver}, counter[, alias])

Defines the addtocounter! function in such a way that a multi-dimensional constraint of type counter is always counted as a single entry. This macro may be called at most once for each counter when defining a solver. As a consequence, extract_sos will may have an Int or a UnitRange{Int} (if multiple constraints are required, which will never be the case for matrix cones) as index parameter for the counter type.

If the alias parameter is present, whenever the type counter is encountered, the counter of type alias is incremented instead. Do not use @counter_alias together with this macro on the same counter.

counter may either be a Symbol, a tuple of Symbols, or the value Any. Note that Any has weakest precedence, irrespective of when the macro is called.

source
PolynomialOptimization.Solver.addtocounter!Function
addtocounter!(state, counters::Counters, ::Val{type}, dim::Integer) where {type} -> Int
addtocounter!(state, counters::Counters, ::Val{type}, num::Integer,
    dim::Integer) where {type} -> UnitRange{Int}

Increments the internal information about how many constraints/variables of a certain type were already added. Usually, a solver implementation does not have to overwrite the default implementation; but it might be useful to do so if some types are for example counted in a single counter or bunch of added elements counts as a single entry. dim is the length of the conic constraint or variable, while num indicates that multiple such constraints or variables of the same type are added. For a list of possible symbols type, see the documentation for extract_sos.

See also Counters, @counter_alias, @counter_atomic.

source

Using this information, an additional implementation may be provided for a faster re-optimization of the same problem:

PolynomialOptimization.Solver.poly_optimizeMethod
poly_optimize(::Val{method}, oldstate, relaxation::AbstractRelaxation,
    groupings::RelaxationGroupings; representation, verbose, kwargs)

A solver that supports re-optimization of an already optimized problem with changed rotations on the DD and SDD representations should implement this method. It is guaranteed that only the rotations change, and only in a structure-preserving way (diagonal and nondiagonal rotations will not interchange). The return value is as documented in poly_optimize, and the oldstate parameter holds the first return value of the previous call to poly_optimize.

source

While this page documents in detail how a new solver can be implemented, the explanation is far more extensive than an actual implementation. In order to implement a new solver, it is therefore recommended to first determine the category in which it falls, then copy and modify an appropriate existing implementation.

AbstractSparseMatrixSolver

This solver type accumulates data in a COO matrix-like format. The callback functions can use append! to quickly add given data to the temporary storage. However, the format is not yet suitable for passing data to the solver, as all monomials are densely indexed. Therefore, in a postprocessing step, the COO matrices have to be converted to CSC matrices with continuous monomial indices using coo_to_csc!. After the optimization is done, the optimal moment vector (the decision variables for moment optimization, the constraint duals for SOS optimization) can be constructed using MomentVector.

PolynomialOptimization.Solver.AbstractSparseMatrixSolverType
AbstractSparseMatrixSolver{I<:Integer,K<:Integer,V<:Real}

Superclass for a solver that requires its data in sparse matrix form. The data is aggregated in COO form using append! and can be converted to CSC form using coo_to_csc!. The type of the indices in final CSC form is I, where the monomials during construction will be represented by numbers of type K. Any type inheriting from this class is supposed to have a field slack::K which is initialized to -one(K) if K is signed or typemax(K) if it is unsigned.

See also SparseMatrixCOO.

source
PolynomialOptimization.Solver.SparseMatrixCOOType
SparseMatrixCOO{I<:Integer,K<:Integer,V<:Real,Offset}

Representation of a sparse matrix in COO form. Fields are rowinds::FastVec{I}, moninds::FastVec{K} (where K is of the type returned by monomial_index), and nzvals::FastVec{V}. The first row/column for the solver has index Offset (of type I).

If this matrix is used in the context of a vectorized COO matrix, the fields are rows, cols, and vals with K === I.

source
Base.append!Method
append!(coo::SparseMatrixCOO, indvals::Indvals)

Appends the data given in indvals into the next row in coo. Returns the index of the row that was added.

See also Indvals.

source
Base.append!Method
append!(coo::SparseMatrixCOO, indvals::IndvalsIterator)

Appends the data given in indvals into successive rows in coo (first(indvals) to the first row, the next to the second, ...). Returns the range of indices of all rows that were added.

See also IndvalsIterator.

source
PolynomialOptimization.Solver.coo_to_csc!Method
coo_to_csc!(coo::Union{SparseMatrixCOO{I,K,V}},Tuple{FastVec{K},FastVec{V}}}...)

Converts sparse COO matrix or vector representations, where the monomial indices of the coo matrices or the entries of the vectors can be arbitrarily sparse, to a CSC-based matrix representation with continuous columns, and the vectors are converted to dense ones. No more than two matrices may be supplied. The input data may be mutated; and this mutated data must be passed on to MomentVector. The following values are returned:

  • number of distinct columns
  • for each input, if it is a matrix, a tuple containing the colptr, rowval, nzval vectors
  • for each input, if it is a vector, the corresponding dense vector
source
PolynomialOptimization.MomentVectorMethod
MomentVector(relaxation::AbstractRelaxation, moments::Vector{<:Real},
    slack::Integer, coo::SparseMatrixCOO...)

Given the moments vector as obtained from a AbstractSparseMatrixSolver solver, convert it to a MomentVector. Note that this function is not fully type-stable, as the result may be based either on a dense or sparse vector depending on the relaxation. To establish the mapping between the solver output and the actual moments, all the column-sorted COO data (i.e., as returned by coo_to_csc!) used in the problem construction needs to be passed on. slack must contain the current value of the slack field of the AbstractSparseMatrixSolver.

source

AbstractAPISolver

This solver type accumulates the data continuously by the means of API calls provided by the solver. The monomials therefore have to have contiguously defined indices from the beginning. This is done by internal bookkeeping; however, it requires the implementation of an additional callback function append! to add new monomials to the solver. After the optimization is done, the optimal moment vector (the decision variables for moment optimization, the constraint duals for SOS optimization) can be constructed using MomentVector.

PolynomialOptimization.Solver.AbstractAPISolverType
AbstractAPISolver{K<:Integer,T,V} <: AbstractSolver{T,V}

Superclass for a solver that requires new variables/constraints to be added via API calls. Solvers that are of this type must implement append! in such a way that they directly add a variable (moment-case) to or constraint (SOS-case) to the solver. Concrete types that inherit from AbstractAPISolver must have a property mon_to_solver::Dict{FastKey{K},T}.

source
Base.append!Method
append!(solver::AbstractAPISolver{K}, key::K)

Appends at least one new variable (moment-case) or constraint (SOS-case) to the solver state that represents the monomial given by key.

source
PolynomialOptimization.MomentVectorMethod
MomentVector(relaxation::AbstractRelaxation, moments::Vector{<:Real},
    solver::AbstractAPISolver)

Given the moments vector as obtained from an AbstractAPISolver, convert it to a MomentVector. Note that this function is not fully type-stable, as the result may be based either on a dense or sparse vector depending on the relaxation.

source

Defining solver capabilities

There are some functions that should be implemented to tell PolynomialOptimization what kind of data the solver expects and which cones are supported; these should return constants.

PolynomialOptimization.Solver.supports_ddFunction
supports_dd(::AbstractSolver)

This function indicates whether the solver natively supports a diagonally-dominant cone (or its dual for the moment case). If it returns false (default), the constraint will be rewritten in terms of multiple $\ell_\infty$/$\ell_1$ norm constraints (if supported, see supports_lnorm), together with slack variables and equality constraints. If these $\ell$-norm constraints are also not supported, linear constraints will be used.

source
PolynomialOptimization.Solver.supports_dd_complexFunction
supports_dd_complex(::AbstractSolver)

This function indicates whether the solver natively supports a complex-valued diagonally-dominant cone (or its dual for the moment case). If it returns false (default), the constraint will be rewritten in terms of quadratic constraints (if supported, see supports_quadratic) or multiple $\ell_\infty$/$\ell_1$ norm constraints (if supported, see supports_lnorm_complex).

source
PolynomialOptimization.Solver.supports_lnormFunction
supports_lnorm(::AbstractSolver)

Indicates the solver support for $\ell_\infty$ (in the moment case) and $\ell_1$ (in the SOS case) norm cones: if true, the cone $x_1 \geq \max_{i \geq 2} \lvert x_i\rvert$ or $x_1 \geq \sum_{i \geq 2} \lvert x_i\rvert$ is supported. The default implementation returns false.

source
PolynomialOptimization.Solver.supports_lnorm_complexFunction
supports_lnorm_complex(::AbstractSolver)

Indicates the solver support for complex-valued $\ell_\infty$ (in the moment case) and $\ell_1$ (in the SOS case) norm cones: if true, the cone $x_1 \geq \max_{i \geq 2} \lvert\operatorname{Re} x_i + \mathrm i \operatorname{Im} x_i\rvert$ or $x_1 \geq \sum_{i \geq 2} \lvert\operatorname{Re} x_i + \mathrm i \operatorname{Im} x_i\rvert$ is supported. The default implementation returns false.

source
PolynomialOptimization.Solver.PSDIndextypeMatrixCartesianType
PSDIndextypeMatrixCartesian(triangle, offset) <: PSDIndextype

The solver implements PSD matrix constraints by using a monolithic PSD matrix variable or an LMI-style representation. Entries from the variable are obtained (or put into the LMI) by using a cartesian index of two integers of the type parameter T of the AbstractSolver. This index represents one triangle of the matrix (the lower if triangle === :L, the upper if triangle === :U). The first entry has the index (offset, offset), typically either 0 or 1.

If this index type is used with primal_moment_setup!, the resulting data will be an iterator through SparseMatrixCOO. In this case, triangle === :F is also permitted, resulting in the full triangle.

Info

Note that even if only one triangle is indexed, it is assumed that the solver will by default populate the other triangle in a completely symmetric way. This corresponds to the typical behavior of solvers that expose PSD variables and allow accessing elements in them via sparse symmetric matrices where only one triangle is given, but the other half is implicit.

See also PSDMatrixCartesian.

source
PolynomialOptimization.Solver.PSDIndextypeVectorType
PSDIndextypeVector(triangle[, scaling]) <: PSDIndextype

The solver implements PSD matrix constraints by demanding that the matrixization of a vector of decision variables be PSD. This index type is not permitted for use with primal_moment_setup!.

If triangle === :F, the vector is formed by stacking all the columns of the matrix. scaling should be omitted.

If triangle === :L, the columns of the lower triangle are assumed to be stacked and scaled, i.e., off-diagonal variables that enter the cone are implicitly multiplied by 1 / scaling in the matrix; so the coefficients will already be premultiplied by scaling (for the add_constr_psd! case) or by 1 / scaling (for the add_var_psd! case). The default value for scaling is $sqrt2$; however, the parameter has to be specified explictly in order to make sure the scaling has the correct type. Note: if no scaling is desired, the preferred value is true, which is equivalent to a multiplicative identity; however, the multiplication can be completely removed during compilation. This is because the value false, which would mean that alll off-diagonal entries are set to zero, is explicitly forbidden.

If triangle === :U, the columns of the upper triangle are assumed to be stacked and scaled.

See also IndvalsIterator.

source
PolynomialOptimization.Solver.PSDIndextypeCOOVectorizedType
PSDIndextypeCOOVectorized(triangle[, scaling], offset)

The solver implements constraints on a monolithic PSD matrix variable using row-by-row scalar products with the vectorized matrix. During vectorization, only the specified triangle is retained. This index type is valid only for the use with primal_moment_setup!.

If triangle === :F, the full matrix is stacked column-wise; scaling should be omitted.

If triangle === :L, the columns of the lower triangled are assumed to be stacked and scaled, i.e., off-diagonal variables that enter the cone are implicitly multiplied by 1 / scaling in the matrix; so the coefficients will already be premultiplied by 1 / scaling. If the solver internally works with the vectorized version, the appropriate value is probably $\sqrt2$; if the solver automatically rewrites everything for full matrices, the appropriate value is either true or $2$.

If triangle === :U, the column of the upper triangle are assumed to be stacked and scaled.

This all refers to the column index of the constraint matrix, where the row index is the index of the constraint. The data is supplied in COO form with specified offset and may be converted to CSC or CSR as desired.

source
PolynomialOptimization.Solver.negate_fixFunction
negate_fix(::AbstractSolver)

Depending on the exact definition of equality constraints (where to put the minus), the dual solutions, i.e., the SOS decomposition, may yield wrong values; then define this function to return true, which will pass the equality constraints with opposite sign.

source
PolynomialOptimization.Solver.negate_freeFunction
negate_free(::AbstractSolver)

Depending on the exact definition of equality constraints (where to put the minus), the dual solutions, i.e., the SOS decomposition, may yield wrong values; then define this function to return true, which will flip the sign of free variables.

source
PolynomialOptimization.Solver.prepend_fixFunction
prepend_fix(::AbstractSolver)

If this method yields true (default), fixed constraints will be created before all others. If it is false, PSD (or nonnegative) constraints will be created first.

Info

Note that this does not imply a global order; if DD or SDD cones are used, fixed constraints may still appear at an arbitrary position. However, this method guarantees the order with respect to cones that will use the same monomial indices.

source
PolynomialOptimization.Solver.prepend_freeFunction
prepend_free(::AbstractSolver)

If this method yields true (default), free variables will be created before all others. If it is false, PSD (or nonnegative) variables will be created first.

Info

Note that this does not imply a global order; if DD or SDD cones are used, free variables may still appear at an arbitrary position. However, this method guarantees the order with respect to cones that will use the same monomial indices.

source

Working with data from the interface

The interface functions that need to be implemented will get their data in the form of index-value pairs or a specific structure related to the desired matrix format.

PolynomialOptimization.Solver.IndvalsType
Indvals{T,V}

Supertype for an iterable that returns a Tuple{T,V} on iteration, where the first is a variable/constraint index and the second its coefficient in the constraint matrix. The parameters of the type correspond to those in AbstractSolver. The properties indices and values can be accessed and will give AbstractVectors of the appropriate type. Note that the fields should only be used if an iterative approach is not feasible, as they might be constructed on-demand (this will only happen for the first two indvals in the standard quadratic cone, all other elements can be accessed with zero cost).

source
PolynomialOptimization.Solver.PSDMatrixCartesianType
PSDMatrixCartesian

An iterable that returns matrix elements of a PSD cone. Iterating through it will give Pairs with keys that contain the mindex of the monomial, and a 3-Tuple of AbstractVectors as the values which contain the row and column indices together with their coefficients to describe where the monomial appears. Note that the vectors will be reused in every iteration.

See also PSDIndextypeMatrixCartesian.

source
PolynomialOptimization.Solver.IndvalsIteratorType
IndvalsIterator{T,V}

An iterable that returns consecutive elements in a vectorized PSD cone. This type stores a vector of indices and values together with information about the length of the individual subsequences. Iterating through it will give Indvals that contain views into the indices and the values. The vector of indices is available via SparseArrays.rowvals, the vector of values via SparseArrays.nonzeros, and the lengths of the subsequences via Base.index_lengths.

See also PSDIndextypeVector.

source

Interface for the moment optimization

The custom implementation of poly_optimize first has to set up all the necessary initial data; then, a call to moment_setup! is sufficient to trigger the process.

PolynomialOptimization.Solver.moment_setup!Function
moment_setup!(state::AbstractSolver, relaxation::AbstractRelaxation,
    groupings::RelaxationGroupings[; representation])

Sets up all the necessary moment matrices, variables, constraints, and objective of a polynomial optimization problem problem according to the values given in grouping (where the first entry corresponds to the basis of the objective, the second of the equality, the third of the inequality, and the fourth of the PSD constraints). The function returns a Vector{<:Vector{<:Tuple{Symbol,Any}}} that contains internal information on the problem. This information is required to obtain dual variables and re-optimize the problem and should be stored in the state.

The following methods must be implemented by a solver to make this function work:

Indices

The variable indices used in all solver functions directly correspond to the indices given back by mindex. However, in a sparse problem there may be far fewer indices present; therefore, when the problem is finally given to the solver, care must be taken to eliminate all unused indices. The functionality provided by AbstractAPISolver and AbstractSparseMatrixSolver already takes care of this.

Order

The individual constraint types can be added in any order (including interleaved).

Representation

This function may also be used to describe simplified cones such as the (scaled) diagonally dominant one. The representation parameter can be used to define a representation that is employed for the individual groupings. This may either be an instance of a RepresentationMethod - which requires the method to be independent of the dimension of the grouping - or a callable. In the latter case, it will be passed as a first parameter an identifier[3] of the current conic variable, and as a second parameter the side dimension of its matrix. The method must then return a RepresentationMethod instance.

See also sos_setup!, moment_add_matrix!, moment_add_equality!, RepresentationMethod.

source
PolynomialOptimization.Solver.moment_add_matrix!Function
moment_add_matrix!(state::AbstractSolver,
    grouping::Union{<:IntMonomialVector,<:AbstractVector{<:IntMonomialVector}},
    constraint::Union{<:IntPolynomial,<:AbstractMatrix{<:IntPolynomial}},
    representation::RepresentationMethod=RepresentationPSD())

Parses a constraint in the moment hierarchy with a basis given in grouping (this might also be a partial basis due to sparsity), premultiplied by constraint (which may be the unit polynomial for the moment matrix) and calls the appropriate solver functions to set up the problem structure according to representation. If constraint is a matrix, grouping should be a vector of monomial vectors (with length corresponding to the side dimension); else, it should just be a monomial vector.

To make this function work for a solver, implement the following low-level primitives:

Usually, this function does not have to be called explicitly; use moment_setup! instead.

See also moment_add_equality!, RepresentationMethod.

source
PolynomialOptimization.Solver.moment_add_equality!Function
moment_add_equality!(state::AbstractSolver, groupings::AbstractVector{<:IntMonomialVector},
    constraint::IntPolynomial)

Parses a polynomial equality constraint for moments and calls the appropriate solver functions to set up the problem structure. groupings contains the bases that will be squared individually in the process to generate the prefactor.

To make this function work for a solver, implement the following low-level primitives:

Usually, this function does not have to be called explicitly; use moment_setup! instead.

See also moment_add_matrix!.

source

For this to work, the following methods (or a subset as previously indicated) must be implemented.

PolynomialOptimization.Solver.add_constr_nonnegative!Function
add_constr_nonnegative!(state::AbstractSolver{T,V}, indvals::Indvals{T,V}) where {T,V}

Add a nonnegative constraint to the solver that contains the decision variables (columns in the linear constraint matrix) indexed according to indvals. Falls back to the vector-valued version if not implemented.

See also Indvals.

source
add_constr_nonnegative!(state::AbstractSolver{T,V},
    indvals::IndvalsIterator{T,V}) where {T,V}

Adds multiple nonnegative constraints to the solver that contain the decision variables (columns in the linear constraint matrix) indices according to the entries in indvals. Falls back to calling the scalar-valued version multiple times if not implemented.

See also IndvalsIterator.

source
PolynomialOptimization.Solver.add_constr_rotated_quadratic!Function
add_constr_rotated_quadratic!(state::AbstractSolver{T,V},
    indvals::IndvalsIterator{T,V}) where {T,V}

Adds a rotated quadratic constraint to the N = length(indvals) linear combinations of decision variables (columns in the conic constraint matrix) indexed according to the indvals. This will read (where $X_i$ is $\mathit{indvals}_i.\mathit{values} \cdot x_{\mathit{indvals}_i.\mathit{indices}}$) $X_1, X_2 \geq 0$, $2X_1 X_2 \geq \sum_{i = 3}^N X_i^2$.

See also Indvals, IndvalsIterator.

Number of parameters

In the real-valued case, indvals is always of length three, in the complex case, it is of length four. If the scaled diagonally dominant representation is requested, indvals can have any length.

Warning

This function will only be called if supports_rotated_quadratic returns true for the given state. If (rotated) quadratic constraints are unsupported, a fallback to a 2x2 PSD constraint is used.

source
PolynomialOptimization.Solver.add_constr_quadratic!Function
add_constr_quadratic!(state::AbstractSolver{T,V}, indvals::IndvalsIterator{T,V}) where {T,V}

Adds a quadratic constraint to the N = length(indvals) linear combinations of decision variables (columns in the conic constraint matrix) indexed according to the indvals. This will read (where $X_i$ is $\mathit{indvals}_i.\mathit{values} \cdot x_{\mathit{indvals}_i.\mathit{indices}}$) $X_1 \geq 0$, $X_1^2 \geq \sum_{i = 2}^N X_i^2$.

See also Indvals, IndvalsIterator.

Number of parameters

In the real-valued case, indvals is always of length three, in the complex case, it is of length four. If the scaled diagonally dominant representation is requested, indvals can have any length.

Warning

This function will only be called if supports_quadratic returns true for the given state. If (rotated) quadratic constraints are unsupported, a fallback to a 2x2 PSD constraint is used.

source
PolynomialOptimization.Solver.add_constr_psd!Function
add_constr_psd!(state::AbstractSolver{T,V}, dim::Integer,
    data::PSDMatrixCartesian{T,V}) where {T,V}

Add a PSD constraint of side dimension dim ≥ 3 to the solver. Its requested triangle is indexed according to the return value of psd_indextype); these elements make up a linear matrix inequality with variables given by the keys when iterating through data, which are of the type T. Note that if add_constr_quadratic! is not implemented, dim may also be 2. This method is called if psd_indextype returns a PSDIndextypeMatrixCartesian.

Complex-valued PSD variables

Note that this function will also be called for complex-valued PSD cones if supports_psd_complex returns false. The data will have been rewritten in terms of a real-valued PSD cone, which doubles the dimension. If the solver natively supports complex-valued PSD cones, add_constr_psd_complex! must be implemented.

source
add_constr_psd!(state::AbstractSolver{T,V}, dim::Integer, data::IndvalsIterator{T,V}) where {T,V}

Add a PSD constraint of side dimension dim ≥ 3 to the solver. data is an iterable through the elements of the PSD matrix one-by-one, in the order specified by psd_indextype. The individual entries are Indvals. This method is called if psd_indextype returns a PSDIndextypeVector.

Complex-valued PSD variables

Note that this function will also be called for complex-valued PSD cones if supports_psd_complex returns false. The data will have been rewritten in terms of a real-valued PSD cone, which doubles the dimension. If the solver natively supports complex-valued PSD cones, add_constr_psd_complex! must be implemented.

source
PolynomialOptimization.Solver.add_constr_psd_complex!Function
add_constr_psd_complex!(state::AbstractSolver{T,V}, dim::Int,
    data::PSDMatrixCartesian{T,Complex{V}}) where {T,V}

Add a Hermitian PSD constraint of side dimension dim ≥ 3 to the solver. Its requested triangle is indexed according to the return value of psd_indextype); these elements make up a linear matrix inequality with variables given by the keys when iterating through data, which are of the type T. The real part of any coefficient corresponds to the coefficient in front of the real part of the matrix entry, the imaginary part is the coefficient for the imaginary part of the matrix entry. Note that if add_constr_quadratic! is not implemented, dim may also be 2. This method is called if psd_indextype returns a PSDIndextypeMatrixCartesian.

Warning

This function will only be called if supports_psd_complex is defined to return true for the given state.

source
add_constr_psd_complex!(state::AbstractSolver{T,V}, dim::Int, data::IndvalsIterator{T,V}) where {T,V}

Add a Hermitian PSD constraint of side dimension dim ≥ 3 to the solver. data is an iterable through the elements of the PSD matrix one-by-one, in the order specified by psd_indextype. The individual entries are Indvals. This method is called if psd_indextype returns a PSDIndextypeVector. Regardless of the travelling order, for diagonal elements, there will be exactly one entry, which is the real part. For off-diagonal elements, the real part will be followed by the imaginary part. Therefore, the coefficients are real-valued.

Warning

This function will only be called if supports_psd_complex is defined to return true for the given state.

source
PolynomialOptimization.Solver.add_constr_dddual!Function
add_constr_dddual!(state::AbstractSolver{T,V}, dim::Integer, data::IndvalsIterator{T,V},
    u) where {T,V}

Add a constraint for membership in the dual cone to diagonally dominant matrices to the solver. data is an iterator through the scaled lower triangle of the matrix. A basis change is induced by u, with the meaning for the primal cone that M ∈ DD(u) ⇔ M = uᵀ Q u with Q ∈ DD.

Warning

This function will only be called if supports_dd returns true for the given state. If diagonally dominant cones are not supported directly, a fallback to a columnwise representation in terms of $\ell_\infty$ norms will be used (or the fallbacks if this norm is not supported).

source
PolynomialOptimization.Solver.add_constr_dddual_complex!Function
add_constr_dddual_complex!(state::AbstractSolver{T,V}, dim::Integer,
    data::IndvalsIterator{T,V}, u) where {T,V}

Add a constraint for membership in the dual cone to complex-valued diagonally dominant matrices to the solver. data is an iterator through the scaled lower triangle of the matrix. A basis change is induced by u, with the meaning for the primal cone that M ∈ DD(u) ⇔ M = u† Q u with Q ∈ DD. For diagonal elements, there will be exactly one entry, which is the real part. For off-diagonal elements, the real part will be followed by the imaginary part. Therefore, the coefficients are real-valued.

Warning

This function will only be called if supports_dd_complex returns true for the given state. If complex-valued diagonally dominant cones are not supported directly, a fallback to quadratic constraints on the complex-valued data is tried first (if supported), followed by a columnwise representation in terms of $\ell_\infty$ norms or their fallback on the realification of the matrix data if not.

source
PolynomialOptimization.Solver.add_constr_linf!Function
add_constr_linf!(state::AbstractSolver{T,V}, indvals::IndvalsIterator{T,V}) where {T,V}

Adds an $\ell_\infty$ norm constraint to the N = length(indvals) linear combinations of decision variables (columns in the conic constraint matrix) indexed according to the indvals. This will read (where $X_i$ is $\mathit{indvals}_i.\mathit{values} \cdot x_{\mathit{indvals}_i.\mathit{indices}}$) $X_1 \geq \max_{i > 2} \lvert X_i\rvert$.

See also Indvals, IndvalsIterator.

Warning

This function will only be called if supports_lnorm returns true for the given state. If $\ell_\infty$ norm constraints are unsupported, a fallback to multiple linear constraints will be used.

source
PolynomialOptimization.Solver.add_constr_linf_complex!Function
add_constr_linf_complex!(state, indvals::IndvalsIterator{T,V}) where {T,V<:Real}

Same as add_constr_linf!, but now two successive items in indvals (starting from the second) are interpreted as determining the real and imaginary part of a component of the $\ell_\infty$ norm cone.

Warning

This function will only be called if supports_lnorm_complex returns true for the given state. If complex-valued $\ell_\infty$ norm constraints are unsupported, a fallback to multiple linear constraints and quadratic cones will be used. If supports_quadratic is not true, complex-valued DD cones cannot be used.

source
PolynomialOptimization.Solver.add_constr_sdddual!Function
add_constr_sdddual!(state::AbstractSolver{T,V}, dim::Integer, data::IndvalsIterator{T,V},
    u) where {T,V}

Add a constraint for membership in the dual cone to scaled diagonally dominant matrices to the solver. data is an iterator through the (unscaled) lower triangle of the matrix. A basis change is induced by u, with the meaning for the primal cone that M ∈ SDD(u) ⇔ M = uᵀ Q u with Q ∈ SDD.

Warning

This function will only be called if supports_sdd returns true for the given state. If scaled diagonally dominant cones are not supported directly, a fallback to (rotated) quadratic cones will be used.

source
PolynomialOptimization.Solver.add_constr_sdddual_complex!Function
add_constr_sdddual_complex!(state::AbstractSolver{T,V}, dim::Integer,
    data::IndvalsIterator{T,V}, u) where {T,V}

Add a constraint for membership in the dual cone to complex-valued scaled diagonally dominant matrices to the solver. data is an iterator through the (unscaled) lower triangle of the matrix. A basis change is induced by u, with the meaning for the primal cone that M ∈ SDD(u) ⇔ M = u† Q u with Q ∈ SDD. For diagonal elements, there will be exactly one entry, which is the real part. For off-diagonal elements, the real part will be followed by the imaginary part. Therefore, the coefficients are real-valued.

Warning

This function will only be called if supports_sdd_complex returns true for the given state. If complex-valued sclaed diagonally dominant cones are not supported directly, a fallback to quadratic constraints is automatically performed.

source
PolynomialOptimization.Solver.add_constr_fix!Function
add_constr_fix!(state::AbstractSolver{T,V}, constrstate, indvals::Indvals{T,V},
    rhs::V) where {T,V}

Add a constraint fixed to rhs to the solver that is composed of all variables (columns in the linear constraint matrix) indexed according to indvals. The parameter constrstate is, upon first call, the value returned by add_constr_fix_prepare!; and on all further calls, it will be the return value of the previous call. Note that rhs will almost always be zero, so if the right-hand side is represented by a sparse vector, it is worth checking for this value (the compiler will be able to remove the check).

See also Indvals.

source
PolynomialOptimization.Solver.add_var_slack!Function
add_var_slack!(state::AbstractSolver{T}, num::Int)

Creates num slack variables in the problem. Slack variables must be free. The result should be an abstract vector (typically a unit range) that contains the indices of all created slack variables. The indices should be of type T.

source

Interface for the SOS optimization

This is in complete analogy to the moment case; the entry point is now sos_setup!.

PolynomialOptimization.Solver.sos_setup!Function
sos_setup!(state::AbstractSolver, relaxation::AbstractRelaxation,
    groupings::RelaxationGroupings[; representation])

Sets up all the necessary SOS matrices, free variables, objective, and constraints of a polynomial optimization problem problem according to the values given in grouping (where the first entry corresponds to the basis of the objective, the second of the equality, the third of the inequality, and the fourth of the PSD constraints). The function returns a Vector{<:Vector{<:Tuple{Symbol,Any}}} that contains internal information on the problem. This information is required to obtain dual constraints and re-optimize the problem and should be stored in the state.

The following methods must be implemented by a solver to make this function work:

Indices

The constraint indices used in all solver functions directly correspond to the indices given back by mindex. However, in a sparse problem there may be far fewer indices present; therefore, when the problem is finally given to the solver, care must be taken to eliminate all unused indices. The functionality provided by AbstractAPISolver and AbstractSparseMatrixSolver already takes care of this.

Order

The individual variable types can be added in any order (including interleaved).

Representation

This function may also be used to describe simplified cones such as the (scaled) diagonally dominant one. The representation parameter can be used to define a representation that is employed for the individual groupings. This may either be an instance of a RepresentationMethod - which requires the method to be independent of the dimension of the grouping - or a callable. In the latter case, it will be passed as a first parameter an identifier[3] of the current conic variable, and as a second parameter the side dimension of its matrix. The method must then return a RepresentationMethod instance.

See also moment_setup!, sos_add_matrix!, sos_add_equality!, RepresentationMethod.

source
PolynomialOptimization.Solver.sos_add_matrix!Function
sos_add_matrix!(state::AbstractSolver, grouping::IntMonomialVector,
    constraint::Union{<:IntPolynomial,<:AbstractMatrix{<:IntPolynomial}},
    representation::RepresentationMethod=RepresentationPSD())

Parses a SOS constraint with a basis given in grouping (this might also be a partial basis due to sparsity), premultiplied by constraint (which may be the unit polynomial for the SOS cone membership) and calls the appropriate solver functions to set up the problem structure according to representation.

To make this function work for a solver, implement the following low-level primitives:

Usually, this function does not have to be called explicitly; use sos_setup! instead.

See also sos_add_equality!.

source
PolynomialOptimization.Solver.sos_add_equality!Function
sos_add_equality!(state::AbstractSolver, groupings::AbstractVector{<:IntMonomialVector},
    constraint::IntPolynomial)

Parses a polynomial equality constraint for sums-of-squares and calls the appropriate solver functions to set up the problem structure. groupings contains the bases that will be squared individually in the process to generate the prefactor.

To make this function work for a solver, implement the following low-level primitives:

Usually, this function does not have to be called explicitly; use sos_setup! instead.

See also sos_add_matrix!.

source

The following methods (or a subset as previously indicated) must be implemented.

PolynomialOptimization.Solver.add_var_nonnegative!Method
add_var_nonnegative!(state::AbstractSolver{T,V}, indvals::Indvals{T,V}) where {T,V}

Add a nonnegative decision variable to the solver and put its value into the linear constraints (rows in the linear constraint matrix) indexed according to indvals. Falls back to the vector-valued version if not implemented.

See also Indvals.

source
PolynomialOptimization.Solver.add_var_nonnegative!Method
add_var_nonnegative!(state::AbstractSolver{T,V}, indvals::IndvalsIterator{T,V}) where {T,V}

Add multiple nonnegative decision variables to the solver and put their values into the linear constraints (rows in the linear constraint matrix) indexed according to the entries in indvals. Falls back to calling the scalar-valued version multiple times if not implemented.

See also IndvalsIterator.

source
PolynomialOptimization.Solver.add_var_rotated_quadratic!Function
add_var_rotated_quadratic!(state::AbstractSolver{T,V},
    indvals::IndvalsIterator{T,V}) where {T,V}

Adds decision variables in a rotated quadratic cone to the solver and put their values into the linear constraints (rows in the linear constraint matrix), indexed according to indvals. The N = length(indvals) variables will satisfy $x_1, x_2 \geq 0$, $2x_1 x_2 \geq \sum_{i = 3}^N x_i^2$.

See also Indvals, IndvalsIterator.

Number of parameters

In the real-valued case, indvals is always of length three, in the complex case, it is of length four. If the scaled diagonally dominant representation is requested, indvals can have any length.

Warning

This function will only be called if supports_quadratic returns true for the given state. If (rotated) quadratic constraints are unsupported, a fallback to a 2x2 PSD variable is used.

source
PolynomialOptimization.Solver.add_var_quadratic!Function
add_var_quadratic!(state::AbstractSolver{T,V}, indvals::IndvalsIterator{T,V}) where {T,V}

Adds decision variables in a quadratic cone to the solver and put their values into the linear constraints (rows in the linear constraint matrix), indexed according to indvals. The N = length(indvals) variables will satisfy $x_1 \geq 0$, $x_1^2 \geq \sum_{i = 2}^N x_i^2$.

See also Indvals, IndvalsIterator.

Number of parameters

In the real-valued case, indvals is always of length three, in the complex case, it is of length four. If the scaled diagonally dominant representation is requested, indvals can have any length.

Warning

This function will only be called if supports_quadratic returns true for the given state. If (rotated) quadratic constraints are unsupported, a fallback to a 2x2 PSD variable is used.

source
PolynomialOptimization.Solver.add_var_psd!Method
add_var_psd!(state::AbstractSolver{T,V}, dim::Int,
    data::PSDMatrixCartesian{T,V}) where {T,V}

Add a PSD variable of side dimension dim ≥ 3 to the solver. Its requested triangle is indexed according to the return value of psd_indextype); these elements of the matrix are put into the linear constraints (rows in the linear constraint matrix) indicated by the keys when iterating through data, which are of the type T, at positions and with coefficients given by their values. Note that if add_var_quadratic! is not implemented, dim may also be 2. This method is called if psd_indextype returns a PSDIndextypeMatrixCartesian.

Complex-valued PSD variables

Note that this function will also be called for complex-valued PSD cones if supports_psd_complex returns false. The data will have been rewritten in terms of a real-valued PSD cone, which doubles the dimension. If the solver natively supports complex-valued PSD cones, add_var_psd_complex! must be implemented.

source
PolynomialOptimization.Solver.add_var_psd!Method
add_var_psd!(state::AbstractSolver{T,V}, dim::Int, data::IndvalsIterator{T,V}) where {T,V}

Conceptually the same as above; but now, data is an iterable through the elements of the PSD variable one-by-one. The individual entries are Indvals. This method is called if psd_indextype returns a PSDIndextypeVector.

Complex-valued PSD variables

Note that this function will also be called for complex-valued PSD cones if supports_psd_complex returns false. The data will have been rewritten in terms of a real-valued PSD cone, which doubles the dimension. If the solver natively supports complex-valued PSD cones, add_var_psd_complex! must be implemented.

source
PolynomialOptimization.Solver.add_var_psd_complex!Function
add_var_psd_complex!(state::AbstractSolver{T,V}, dim::Int, data::PSDMatrixCartesian{T,Complex{V}}) where {T,V}

Add a Hermitian PSD variable of side dimension dim ≥ 3 to the solver. Its requested triangle is indexed according to the return value of psd_indextype); these elements of the matrix are put into the linear constraints (rows in the linear constraint matrix) indicated by the keys when iterating through data, which are of the type T, at positions and with coefficients given by their values. The real part of the coefficient corresponds to the coefficient in front of the real part of the matrix entry, the imaginary part is the coefficient for the imaginary part of the matrix entry. Note that if add_var_quadratic! is not implemented, dim may also be 2. This method is called if psd_indextype returns a PSDIndextypeMatrixCartesian.

Warning

This function will only be called if supports_psd_complex is defined to return true for the given state.

source
add_var_psd_complex!(state::AbstractSolver{T,V}, dim::Int,
    data::IndvalsIterator{T,V}) where {T,V}

Conceptually the same as above; but now, data is an iterable through the elements of the PSD variable one-by-one. The individual entries are Indvals. This method is called if psd_indextype returns a PSDIndextypeVector. Regardless of the travelling order, for diagonal elements, there will be exactly one entry, which is the real part. For off-diagonal elements, the real part will be followed by the imaginary part. Therefore, the coefficients are real-valued.

Warning

This function will only be called if supports_psd_complex is defined to return true for the given state.

source
PolynomialOptimization.Solver.add_var_dd!Function
add_var_dd!(state::AbstractSolver{T,V}, dim::Integer, data::IndvalsIterator{T,V},
    u) where {T,V}

Add a constraint for membership in the cone of diagonally dominant matrices to the solver. data is an iterator through the scaled lower triangle of the matrix. A basis change is induced by u, with the meaning that M ∈ DD(u) ⇔ M = uᵀ Q u with Q ∈ DD.

Warning

This function will only be called if supports_dd returns true for the given state. If diagonally dominant cones are not supported directly, a fallback to a columnwise representation in terms of $\ell_1$ norms will be used (or the fallbacks if this norm is not supported).

source
PolynomialOptimization.Solver.add_var_dd_complex!Function
add_var_dd_complex!(state::AbstractSolver{T,V}, dim::Integer,
    data::IndvalsIterator{T,V}, u) where {T,V}

Add a constraint for membership in the cone of complex-valued diagonally dominant matrices to the solver. data is an iterator hrough the scaled lower triangle of the matrix. A basis change is induced by u, with the meaning that M ∈ DD(u) ⇔ M = u† Q u with Q ∈ DD. For diagonal elements, there will be exactly one entry, which is the real part. For off-diagonal elements, the real part will be followed by the imaginary part. Therefore, the coefficients are real-valued.

Warning

This function will only be called if supports_dd_complex returns true for the given state. If complex-valued diagonally dominant cones are not supported directly, a fallback to quadratic cones on the complex-valued data is tried first (if supported), followed by a columnwise representation in terms of $\ell_1$ norms or their fallback on the realification of the matrix data if not.

source
PolynomialOptimization.Solver.add_var_l1!Function
add_var_l1!(state::AbstractSolver{T,V}, indvals::IndvalsIterator{T,V}) where {T,V}

Adds decision variables in an $\ell_1$ norm cone to the solver and put their values into the linear constraints (rows in the linear constraint matrix), indexed according to the indvals. The N = length(indvals) variables will satisfy $x_1 \geq \sum_{i = 2}^N \lvert x_i\rvert$.

See also Indvals, IndvalsIterator.

Warning

This function will only be called if supports_lnorm returns true for the given state. If $\ell_\infty$ norm cones are unsupported, a fallback to multiple nonnegative variables will be used.

source
PolynomialOptimization.Solver.add_var_l1_complex!Function
add_var_l1_complex!(state::AbstractSolver{T,V}, indvals::IndvalsIterator{T,V}) where {T,V}

Same as add_var_l1!, but now two successive items in indvals (starting from the second) are interpreted as determining the real and imaginary part of a component of the $\ell_1$ norm variable.

Warning

This function will only be called if supports_lnorm_complex returns true for the given state. If complex-valued $\ell_1$ norm cones are unsupported, a fallback to multiple nonnegative and quadratic variables will be used.

source
PolynomialOptimization.Solver.add_var_sdd!Function
add_var_sdd!(state::AbstractSolver{T,V}, dim::Integer, data::IndvalsIterator{T,V},
    u) where {T,V}

Add a constraint for membership in the cone of scaled diagonally dominant matrices to the solver. data is an iterator through the (unscaled) lower triangle of the matrix. A basis change is induced by u, with the meaning that M ∈ SDD(u) ⇔ M = uᵀ Q u with Q ∈ SDD.

Warning

This function will only be called if supports_sdd returns true for the given state. If scaled diagonally dominant cones are not supported directly, a fallback to (rotated) quadratic cones will be used.

source
PolynomialOptimization.Solver.add_var_sdd_complex!Function
add_var_sdd_complex!(state::AbstractSolver{T,V}, dim::Integer,
    data::IndvalsIterator{T,V}, u) where {T,V}

Add a constraint for membership in the cone of complex-valued scaled diagonally dominant matrices to the solver. data is an iterator through the (unscaled) lower triangle of the matrix. A basis change is induced by u, with the meaning that M ∈ SDD(u) ⇔ M = u† Q u with Q ∈ SDD. For diagonal elements, there will be exactly one entry, which is the real part. For off-diagonal elements, the real part will be followed by the imaginary part. Therefore, the coefficients are real-valued.

Warning

This function will only be called if supports_sdd_complex returns true for the given state. If complex-valued scaled diagonally dominant cones are not supported directly, a fallback to quadratic cones is automatically performed.

source
PolynomialOptimization.Solver.add_var_free!Function
add_var_free!(state::AbstractSolver{T,V}, eqstate, indvals::Indvals{T,V},
    obj::V) where {T,V}

Add a free variable to the solver and put its value into the linear constraints (rows in the linear constraint matrix), indexed according to indvals. The variable should also be put into the objective with coefficient obj (which is likely to be zero). The parameter eqstate is, upon first call, the value returned by add_var_free_prepare!; and on all further calls, it will be the return value of the previous call.

See also Indvals.

source
PolynomialOptimization.Solver.fix_constraints!Function
fix_constraints!(state::AbstractSolver{T,V}, indvals::Indvals{T,V}) where {T,V}

Ensures that all constraints in the optimization problem are fixed to the values according to indvals. This function will be called exactly once by sos_setup! after all variables and constraints have been set up.

See also Indvals.

source
fix_constraints(state::AbstractSolver{<:Integer,V}, m::Int,
    indvals::Indvals{I,V}) where {I,V}

Ensures that all constraints in the optimization problem are fixed to the values according to indvals. This form of the function is called from primal_moment_setup!. m is the number of constraints in the solver.

source
PolynomialOptimization.Solver.add_constr_slack!Function
add_constr_slack!(state::AbstractSolver{T}, num::Int)

Creates num linear fix-to-zero slack constraints in the problem (i.e., constraints that do not correspond to moments). The result should be an abstract vector (typically a unit range) that contains the indices of type T of all created slack constraints.

source

Interface for the moment optimization in primal form

A very particular case is if the moment optimization should be done using semidefinite variables that the solver allows to define, but only as a single variable (not as a cone in which to put variables); this is the primal form, most suitable for SOS optimizations. However, there can be a good reason to use the primal form for moment optimizations instead: namely, if the solver can exploit a low-rank assumptions on this matrix. In this case, poly_optimize should call primal_moment_setup! instead:

PolynomialOptimization.Solver.primal_moment_setup!Function
primal_moment_setup!(state::AbstractSolver, relaxation::AbstractRelaxation,
    groupings::RelaxationGroupings; verbose=false)

Sets up all the necessary moment matrices, constraints, and objective of a polynomial optimization problem in primal form, i.e., for solvers which allow to declare monolithic semidefinite and nonnegative variables and linear constraints. While usually, the SOS form would be more suitable for such solvers, forcing the moment form in primal representation ensures that low-rank assumptions about the primal variable hold true, which can be exploited by some solvers. This function returns a Vector{<:Vector{<:Tuple{Symbol,Any}}} that contains internal information on the problem. This information is required to obtain dual variables and re-optimize the problem and should be stored in the state. It also returns an internal state that is important for the reconstruction of the moment matrices using MomentVector. The internal state has three properties that are part of the public interface:

  • num_con::Int is the number of constraints
  • num_nonneg::Int is the number of nonnegative variables
  • psd_dim::FastVec{I} holds the side dimensions of the PSD variables (where I is the index type of state)

The following methods must be implemented by a solver to make this function work:

Warning

During the reformulation, the function is able to detect a certain class of unbounded or infeasible problem formulations. If this is the case, it will return missing without invoking the solver.

source

The following methods must be implemented:

PolynomialOptimization.Solver.add_var_nonnegative!Method
add_var_nonnegative!(state::AbstractSolver{<:Integer,V}, m::Int, n::Int,
    data::SparseMatrixCOO{I,I,V}, obj::Tuple{FastVec{I},FastVec{V}}) where {I,V}

This form of the function is called from primal_moment_setup!. m is the number of total constraints, n the number of nonnegative variables; and the linear constraint matrix is given by data. The variables should be put into the objective according to obj. Note that the indices in data will have an offset as defined by the psd_indextype, while the indices in obj take their offset from objective_indextype.

See also coo_to_csc!

source
PolynomialOptimization.Solver.add_var_psd!Method
add_var_psd!(state::AbstractSolver{<:Integer,V}, m::Int, dim::I,
    data::SparseMatrixCOO{I,I,V},
    obj::Union{Nothing,Tuple{FastVec{I},FastVec{V}}}) where {I,V}

This form of the function is called from primal_moment_setup! when both the psd_indextype and the objective_indextype are of type PSDIndextypeCOOVectorized. Note that the triangles or scaling factors are allowed to be different. m is the number of constraints in the problem, and the first index in data always corresponds to the constraint index. The matrix also takes part in the objective unless the obj parameter is nothing.

There are various possible variations in which the data can be passed, as documented for the following functions.

See also coo_to_csc!

source
PolynomialOptimization.Solver.fix_constraints!Method
fix_constraints(state::AbstractSolver{<:Integer,V}, m::Int,
    indvals::Indvals{I,V}) where {I,V}

Ensures that all constraints in the optimization problem are fixed to the values according to indvals. This form of the function is called from primal_moment_setup!. m is the number of constraints in the solver.

source

Note that these methods work with the COO representation, which can be quickly converted to either CSR or CSC using coo_to_csr! and coo_to_csc!, which respects the offset desired by the solver. Only this interface allows to set psd_indextype to a PSDIndextypeCOOVectorized; but PSDIndextypeVector is now forbidden.

Helper functions

The solver module exports a number of helper functions which may be of use in implementations:

PolynomialOptimization.Solver.count_uniquesFunction
count_uniques(vec::AbstractVector[, callback])
count_uniques(vec₁::AbstractVector, vec₂::AbstractVector[, callback])

Return the number of unique elements in the vector(s), which must be sorted but may possibly contain duplicates. The callback is invoked once for every unique entry. Its first parameter is the index of the element in the unique total vector, its second (and third) is the last index/indices correponding to the element in the input vector(s). In the second form which allows to check for two vectors jointly, one of the callback parameters can be missing if the element is present only in one of the two vectors.

source
PolynomialOptimization.Solver.coo_to_csc!Method
coo_to_csc!(nCols, s::SparseMatrixCOO)

Converts a COO matrix into a CSC matrix, where the three vectors (colptr, rowvals, nzvals) are returned. The offset of s is respected. Note that some of the vectors in s are modified by this function.

source
PolynomialOptimization.Solver.coo_to_csr!Function
coo_to_csr!(nRows, s::SparseMatrixCOO)

Converts a COO matrix into a CSR matrix, where the three vectors (rowptr, colvals, nzvals) are returned. The offset of s is respected. Note that some of the vectors in s are modified by this function.

source

Additionally, Solver reexports a number of useful types and functions for developing the interface (see src/optimization/solver/Solver.jl); therefore, usually only the Solver submodule itself has to be used and not PolynomialOptimization itself. However, note that it is highly recommended to say using PolynomialOptimization: @assert, @inbounds in the solver implementation; this will replace the Base implementations of the macros by ones that, depending on a debugging constant, enable or disable the desired functionality.

Newton.halfpolytope

Finding the Newton halfpolytope requires a linear solver that supports problem modification for quick reoptimization. All the functions here are defined in the submodule PolynomialOptimization.Newton and they are not exported.

List of supported solvers

The following list contains all the solvers and the required packages that provide access to the solver. The name of the solver is identical with the solver method (as a Symbol).

SolverPackageLicenseSpeedAccuracyMemory
COPTCOPT.jlcommercial👍👍👍👍👍👍👍👍👍
MosekMosek.jlcommercial👍👍👍👍👍👍👍👍👍

Solver interface

The following functions need to be implemented so that a solver is available via Newton.halfpolytope. The preprocessing functions can be omitted if no preprocessing capabilities should be provided.

PolynomialOptimization.Newton.preprocFunction
preproc(V, mons::IntMonomialVector, vertexindices, verbose::Bool, singlethread::Bool;
    parameters...) -> AbstractVector{Bool}

Checks for convex dependencies of exponents from the list of monomials. vertexindices might either be a Vector{Int} that indexes mons or is it Val(:all). In the former case, the convex polytope is spanned by the exponents of all monomials in mons[vertexindices] (and it can be assumed that those are independent, and that vertexindices is fairly short compared to mons); in the latter case, the convex polytope is potentially spanned by the exponents of all monomials in mons. The implementation has to return an AbstractVector{Bool} that for every monomial in mons indicates whether the convex hull spanned by the elements just described (excluding the monomial in question) already contains the monomial. An entry must be false if and only if this is possible, i.e., it is redundant. singlethread will be true if this function will be called in parallel by multiple threads, so that the linear solver itself should be single-threaded.

source
PolynomialOptimization.Newton.prepareFunction
prepare(V, vertices::IntMonomialVector, num::Int, verbose::Bool; parameters...) ->
    (nthreads::Int, userdata, userdata_clone)

This function is responsible for creating an optimization task that can be used to check membership in the Newton halfpolytope. The vertices of the polytope are given by the exponents of mons. The total number of monomials that have to be checked is given by num. The function must return the number of threads that will be used to carry out the optimization (which is not the number of threads that the optimizer uses internally, but determines how PolynomialOptimization will distribute the jobs), an internal optimization task that is passed on to work, and a copy of this task for use in a second thread if the number of threads is greater than one, else nothing. More copies will be created as required by clonetask; however, assuming that setting up the task will potentially require more resources than cloning a task allows the function to estimate the required memory (and therefore a sensible number of threads) better by already performing one clone.

See also @allocdiff.

source
PolynomialOptimization.Newton.alloc_globalFunction
alloc_global(V, nv)

This function is called once in the main thread before work is executed (which, due to multithreading, might occur more than once). It can create some shared data that is used in a read-only manner by all workers at the same time. The default implementation of this function does nothing.

source
PolynomialOptimization.Newton.alloc_localFunction
alloc_local(V, nv)

This function is called once in every computation thread before work is executed (which, due to task splitting, might occur more than once). It can create some shared data that is available for reading and writing by every worker. The default implementation of this function does nothing.

source
PolynomialOptimization.Newton.workFunction
work(V, task, data_global, data_local, expiter::IntMonomialVectorIterator{true},
    Δprogress::Ref{Int}, Δacceptance::Ref{Int}, add_callback::Function,
    iteration_callback::Union{Nothing,Function})

Iterates through expiter (which gives a tuple consisting of the index and the exponents) and for every exponent checks whether this it can be obtained by a convex combination of the coefficients as they are set up in task. If yes, add_callback must be called with the exponent index as a parameter, and Δacceptance should be incremented. In any case, Δprogress should be incremented. Additionally, after every check, iteration_callback should be called with no parameters, if it is a function. The data parameters contain the custom data that was previously generated using alloc_global and alloc_local. Only data_local may be mutated.

source

Once a solver has been implemented, it should add its solver symbol to the vector Newton.newton_methods, which enables this solver to be chosen automatically.

Automatic tightening

Automatic tightening of a polynomial optimization problem requires a linear solver that finds a solution to a system of linear equations that minimizes the ℓ₁ norm (better yet, the ℓ₀-norm, if you can implement this). The solver is only called if the number of rows is smaller than the number of columns; else, the solution is calculated using SPQR's direct solver.

List of supported solvers

The following list contains all the solvers and the required packages that provide access to the solver. The name of the solver is identical with the solver method (as a Symbol).

SolverPackageLicenseSpeedAccuracyMemory
COPTCOPT.jlcommercial👍👍👍👍👍👍👍👍👍
MosekMosek.jlcommercial👍👍👍👍👍👍👍👍👍

Solver interface

The following function needs to be implemented so that a solver is available via automatic tightening.

Once a solver has been implemented, it should add its solver symbol to the vector PolynomialOptimization.tightening_methods, which enables this solver to be chosen automatically.

  • 1If the return type is a vector, psd_indextype should be defined on state, and it must return a PSDIndextypeVector. The scaling property of the index type will automatically be inverted, so that what was $\sqrt2$ before now becomes $\frac{1}{\sqrt2}$.
  • 2Complex values can be treated either by returning a vector of Complex element type, or by returning a real-valued vector where the diagonals (PSD/DD/SDD)/first elements ($\ell$-norm) have a single entry and off-diagonals two.
  • 3This identifier will be a tuple, where the first element is a symbol - either :objective, :nonneg, or :psd - to indicate the general reason why the variable is there. The second element is an Int denoting the index of the constraint (and will be undefined for the objective, but still present to avoid extra compilation). The last element is an Int denoting the index of the grouping within the constraint/objective.