API reference

ClusteredLowRankSolver.BlockType
Block(l::Any[, r::Int, s::Int])

Specifies a block corresponding to the positive semidefinite variable l.

Specifying r,s makes the Block correspond to the r,s subblock of the variable l.

ClusteredLowRankSolver.ClusteredLowRankSDPMethod
ClusteredLowRankSDP(problem::Problem; prec=precision(BigFloat), verbose=false)

Define a ClusteredLowRankSDP based on problem.

Keyword arguments:

  • prec (default: precision(BigFloat)): the precision of the result
  • verbose (default: false): print progress to the standard output
ClusteredLowRankSolver.ConstraintType
Constraint(constant, matrixcoeff, freecoeff[, samples, scalings])

Models a polynomial constaint of the form

\[ ∑_i ⟨ A_i(x), Y_i ⟩ + ∑_j b_j(x) y_j = c(x)\]

using sampling on samples. Here the samples are only required if $A_i$, $b_j$, and $c$ are polynomials. When the coefficient matrix $A_i$ has block structure with equal sized blocks, the Block struct can be used as key to indicate to which subblock the given matrix corresponds.

Arguments:

  • constant: The right hand side $c(x)$
  • matrixcoeff::Dict: The coefficient matrices for the positive semidefinite matrix variables.
  • freecoeff::Dict: The coefficients for the free variables.
  • samples::Vector: The sample set on which the constraint should be satisfied. Required for polynomial constraints.
  • scalings::Vector: Optional; scale the constraint with a factor depending on the sample index.
ClusteredLowRankSolver.DualSolutionType
DualSolution{T}

A dual solution to the semidefinite program, with fields

  • base_ring
  • matrixvars::Dict{Any, Matrix{T}}
  • freevars::Dict{Any, T}
ClusteredLowRankSolver.LowRankMatPolType
LowRankMatPol(lambda::Vector, vs::Vector{Vector}[, ws::Vector{Vector}])

The matrix $∑_i λ_i v_i w_i^{\sf T}$ where $v_i$ are the entries of vs and $w_i$ the entries of ws

If ws is not specified, use ws = vs.

ClusteredLowRankSolver.ObjectiveType
Objective(constant, matrixcoeff::Dict, freecoeff::Dict)

The objective for the Problem.

Arguments:

  • constant: A constant offset of the objective value.
  • matrixcoeff: A Dict with positive semidefinite matrix variable names as keys and the objective matrices as values.
  • freecoeff: A Dict with free variable names as keys and the coefficients as values.
ClusteredLowRankSolver.ProblemType
    Problem(maximize::Bool, objective::Objective, constraints::Vector{Constraint})
    Problem(obj::Union{Maximize, Minimize}, constraints::Vector)

Combine the objective and constraints into a low-rank polynomial problem.

ClusteredLowRankSolver.RoundingSettingsType
RoundingSettings(settings...)

Settings for the rounding procedure:

  1. Finding the kernel
    • kernel_errbound: (default: 1e-10) the allowed error for the kernel vectors. That is, the maximum entry of Xv is in absolute value at most this
    • kernel_round_errbound: (default: 1e-15) the maximum allowed error when rounding the reduced row-echelon form entrywise
    • kernel_use_primal: (default: true) use the primal solution to find the kernel vectors (otherwise, use an SVD of the dual solution)
    • kernel_lll: (default: false) if true, use the LLL algorithm to find the nullspace of the kernel. Otherwise, use the reduced row-echelon form to find kernel vectors.
    • kernel_bits: (default: 1000) the maximum number of bits to be used in the LLL algorithm (for finding relations or finding the entries of the RREF)
  2. Reducing the kernel vectors
    • reduce_kernelvectors: (default: true) apply the reduction step or not
    • reduce_kernelvectors_cutoff: (default: 400) do reduction on the full matrix if its size is at most this cutoff. Otherwise do it on a submatrix
    • reduce_kernelvectors_stepsize: (default: 200) the number of extra columns to take into account in each iteration of the reduction step
  3. Transforming the problem and the solution
    • unimodular_transform: (default: true) use a unimodular transform obtained in the reduction step
    • normalize_transformation: (default: true) multiply by a diagonal matrix to get an integer transformation for the problem (for problems over QQ)
  4. Finding an approximate solution in the field
    • regularization: (default: 1e-20) use this regularization for solving the extended system
    • approximation_decimals: (default: 40) Approximate the numerical solution using this many digits, entrywise
  5. Rounding the solution to the affine space of constraints
    • redundancyfactor: (default: 10) take at least this times the number of constraints columns as potential pivots
    • pseudo: (default: true) use the psuedo inverse for rounding (this may give solutions with larger bitsize than needed)
    • pseudo_columnfactor: (default: 1.05) For a system of r rows, use r * pseudo_columnfactor number of columns for the pseudo inverse
    • extracolumns_linindep: (default: false) if true, take the extra columns linearly independent of each other (otherwise, random columns)
ClusteredLowRankSolver.SampledMPolyRingType
SampledMPolyRing(base_ring, samples)

A sampled polynomial ring with evaluations in base_ring, only defined on the samples. The samples should be sorted.

ClusteredLowRankSolver.SampledMPolyRingElemType
SampledMPolyRinElem(parent::SampledMPolyRing, evaluations::Vector)

A sampled polynomial corresponding to the given parent, which evaluates to evaluations on the samples of the parent ring.

ClusteredLowRankSolver.approximatefeketeMethod
approximatefekete(basis, samples; base_ring=BigFloat) -> basis, samples

Compute approximate fekete points based on samples and a corresponding orthogonal basis. The basis consists of sampled polynomials, sampled at samples.

This preserves a degree ordering of basis if present.

ClusteredLowRankSolver.basis_gegenbauerMethod
basis_gegenbauer(d, n, x)

Basis for the Gegenbauer polynomials in dimension n up to degree d. This is the Gegenbauer polynomial with parameter lambda = n/2-1, or the Jacobi polynomial with alpha = beta = (n-3)/2. Normalized to evaluate to 1 at 1. Taken from arxiv/2001.00256, ancillary files, SemidefiniteProgramming.jl.

ClusteredLowRankSolver.check_sdp!Method
check_sdp!(sdp::ClusteredLowRankSDP)

Check whether the constraint matrices are symmetric, and remove empty constraint matrices.

ClusteredLowRankSolver.convert_to_precFunction
convert_to_prec(sdp, prec=precision(BigFloat))

Convert the semidefinite program to a semidefinite program with prec bits of precision, without error bounds.

ClusteredLowRankSolver.exact_solutionMethod
exact_solution(problem::Problem, primalsol::PrimalSolution, dualsol::DualSolution; transformed=false, FF = QQ, g=1, settings::RoundingSettings=RoundingSettings(), monomial_bases=nothing)

Compute and return an exact solution to the problem, given a primal solution, dual solution and a field FF with approximate generator g. Return (success, exactdualsol) if transformed=false, and (success, pd_transformed_exactsolution, transformations) if transformed=true.

ClusteredLowRankSolver.find_fieldMethod
find_field(primalsol, dualsol, max_degree=4; valbound=1e-15, errbound=1e-15, bits=100, max_coeff=1000)

Heuristically find a field over which the kernel can probably be defined.

Only consider values at least valbound in absolute value. Find minimal polynomials such that the chosen entries are approximately generators with an error bound of errbound. Use bits number of bits and reject minimal polynomials with a maximum coefficient of more than max_coeff.

ClusteredLowRankSolver.generic_embeddingMethod
generic_embedding(exactvalue, g; base_ring=BigFloat)

Convert the exact numbers from a number field to floating point approximations, using a floating point approximation of a generator g of the field.

Convert rationals and integers to the same numbers in base_ring.

ClusteredLowRankSolver.linearsystemMethod
linearsystem(problem::Problem, FF=QQ)

Let x be the vcat of the vectorizations of the matrix variables. This function returns the matrix A and vector b such that the constraints are given by the system Ax = b (over the field FF).

ClusteredLowRankSolver.linearsystem_coefficientmatchingMethod
linearsystem_coefficientmatching(problem::Problem, monomial_basis::Vector{Vector{MPolyRingElem}}; FF=QQ)

Let x be the vcat of the vectorizations of the matrix variables. This function returns the matrix A and vector b such that the constraints obtained from coefficient matching are given by the system Ax = b over the field FF, with one constraint per monomial in monomial_basis. The problem should not contain SampledMPolyRingElem's for this to work.

ClusteredLowRankSolver.model_psd_variables_as_free_variablesMethod
model_psd_variables_as_free_variables!(problem::Problem, as_free)

Model the positive semidefinite variables with names in as_free using free variables, and add extra constraints to set them equal to auxilliary positive semidefinite variables.

ClusteredLowRankSolver.objvalueMethod
objvalue(problem::Problem, sol::DualSolution)

Return the objective value of the dual solution with respect to the given objective or problem.

ClusteredLowRankSolver.partial_linearsystemMethod
partial_linearsystem(problem::Problem, sol::DualSolution, columns::Union{DualSolution, Vector{Int}})

Let x be the vcat of the vectorizations of the matrix variables. Let I be the index set of the columns. This function returns the matrix AI and vector b-Ax such that the constraints for the error vector e using variables indexed by I are given by the system AIe = b-Ax.

ClusteredLowRankSolver.sdpa_sparse_to_problemFunction
sdpa_sparse_to_problem(filename,obj_shift = 0; T=Float64)

Define the Problem from the file filename assuming it is in SDPA sparse format, using the number type T. Optionally add an objective shift.

ClusteredLowRankSolver.slacksMethod
slacks(problem, dualsol)

Compute the difference between the left hand side and the right hand side of all constraints

ClusteredLowRankSolver.solvesdpMethod
	solvesdp(problem::Problem; kwargs...)
    solvesdp(sdp::ClusteredLowRankSDP; kwargs...)

Solve the semidefinite program generated from problem or sdp.

Keyword arguments:

  • prec (default: precision(BigFloat)): the precision used
  • gamma (default: 0.9): the step length reduction; a maximum step length of α reduces to a step length of max(gamma*α,1)
  • beta_(in)feasible (default: 0.1 (0.3)): the amount mu is tried to be reduced by in each iteration, for (in)feasible solutions
  • omega_p/d (default: 10^10): the starting matrix variable for the primal/dual is omega_p/d*I
  • maxiterations (default: 500): the maximum number of iterations
  • duality_gap_threshold (default: 10^-15): how near to optimal the solution needs to be
  • primal/dual_error_threshold (default:10^-30): how feasible the primal/dual solution needs to be
  • max_complementary_gap (default: 10^100): the maximum of dot(X,Y)/nrows(X) allowed
  • need_primal_feasible/need_dual_feasible (default: false): terminate when the solution is primal/dual feasible
  • verbose (default: true): print information after every iteration if true
  • step_length_threshold (default: 10^-7): the minimum step length allowed
  • primalsol (default: nothing): start from the solution (primalsol, dualsol) if both primalsol and dualsol are given
  • dualsol (default: nothing): start from the solution (primalsol, dualsol) if both primalsol and dualsol are given
ClusteredLowRankSolver.to_fieldMethod
to_field(v, N, g; bits=100, errbound=1e-15)

Find an approximation of v in the number field N, using the approximate generator g of N.

ClusteredLowRankSolver.vectorizeMethod
vectorize(sol)

Vectorize the solution by taking the upper triangular part of the matrices. The variables are first sorted by size and then by hash.