Complex Polynomial Optimization
A general complex polynomial optimization problem could be formulized as
\[\mathrm{inf}_{\mathbf{z}\in\mathbf{K}}\ f(\mathbf{z},\bar{\mathbf{z}}),\]
with
\[\mathbf{K}\coloneqq\lbrace \mathbf{z}\in\mathbb{C}^n \mid g_i(\mathbf{z},\bar{\mathbf{z}})\ge0, i=1,\ldots,m,\ h_j(\mathbf{z},\bar{\mathbf{z}})=0, j=1,\ldots,\ell\rbrace,\]
where $\bar{\mathbf{z}}$ stands for the conjugate of $\mathbf{z}:=(z_1,\ldots,z_n)$, and $f, g_i, i=1,\ldots,m, h_j, j=1,\ldots,\ell$ are real-valued complex polynomials satisfying $\bar{f}=f$ and $\bar{g}_i=g_i$, $\bar{h}_j=h_j$.
Let us consider the following example:
\[\mathrm{inf}\ 3-|z_1|^2-0.5\mathbf{i}z_1\bar{z}_2^2+0.5\mathbf{i}z_2^2\bar{z}_1\]
\[\mathrm{s.t.}\ z_2+\bar{z}_2\ge0,|z_1|^2-0.25z_1^2-0.25\bar{z}_1^2=1,|z_1|^2+|z_2|^2=3,\mathbf{i}z_2-\mathbf{i}\bar{z}_2=0.\]
using DynamicPolynomials
n = 2 # set the number of complex variables
@complex_polyvar z[1:n]
f = 3 - z[1]*conj(z[1]) - 0.5im*z[1]*conj(z[2])^2 + 0.5im*z[2]^2*conj(z[1])
g1 = z[2] + conj(z[2])
g2 = z[1]*conj(z[1]) - 0.25*z[1]^2 - 0.25*conj(z[1])^2 - 1
g3 = z[1]*conj(z[1]) + z[2]*conj(z[2]) - 3
g4 = im*z[2] - im*conj(z[2])
pop = [f, g1, g2, g3, g4]
order = 2 # set the relaxation order
opt,sol,data = complex_tssos_first(pop, z, order, numeq=3, TS="block", solution=true) # no correlative sparsity
opt,sol,data = complex_cs_tssos_first(pop, z, order, numeq=3, TS="block", solution=true)
Keyword arguments
Argument | Description | Default value |
---|---|---|
nb | Specify the first nb complex variables to be of unit norm | 0 |
numeq | Specify the last numeq constraints to be equality constraints | 0 |
CS | Types of chordal extensions in exploiting correlative sparsity: "MF" (approximately smallest chordal extension), "NC" (not performing chordal extension), false (invalidating correlative sparsity exploitation) | "MF" |
cliques | Use customized variable cliques | [] |
TS | Types of chordal extensions used in term sparsity iterations: "block"(maximal chordal extension), "MD" (approximately smallest chordal extension), false (invalidating term sparsity iterations) | "block" |
ConjugateBasis | include conjugate variables in monomial bases | false |
normality | Impose normality condtions of order normality | 1 |
merge | Merge overlapping PSD blocks | false |
md | Parameter for tunning the merging strength | 3 |
MomentOne | add a first-order moment PSD constraint for each variable clique | false |
solver | Specify an SDP solver: "Mosek" or "COSMO" | "Mosek" |
cosmo_setting | Parameters for the COSMO solver: cosmo_para(eps_abs, eps_rel, max_iter, time_limit) | cosmo_para(1e-5, 1e-5, 1e4, 0) |
mosek_setting | Parameters for the Mosek solver: mosek_para(tol_pfeas, tol_dfeas, tol_relgap, time_limit, num_threads) | mosek_para(1e-8, 1e-8, 1e-8, -1, 0) |
QUIET | Silence the output | false |
solve | Solve the SDP relaxation | true |
dualize | Solve the dual SDP problem | false |
Gram | Output Gram matrices | false |
solution | Extract an optimal solution | false |
rtol | tolerance for rank | 1e-2 |
gtol | tolerance for global optimality gap | 1e-2 |
ftol | tolerance for feasibility | 1e-3 |
References
- Exploiting Sparsity in Complex Polynomial Optimization, Jie Wang and Victor Magron, 2021.
Methods
Missing docstring for complex_tssos_first
. Check Documenter's build log for details.
Missing docstring for complex_tssos_higher!
. Check Documenter's build log for details.
TSSOS.complex_cs_tssos_first
— Functionopt,sol,data = complex_cs_tssos_first(pop, z, d; nb=0, numeq=0, CS="MF", cliques=[], TS="block", reducebasis=false,
merge=false, md=3, solver="Mosek", QUIET=false, solve=true, solution=false, dualize=false, Gram=false,
MomentOne=false, ConjugateBasis=false, normality=!ConjugateBasis, cosmo_setting=cosmo_para(), mosek_setting=mosek_para(),
rtol=1e-2, gtol=1e-2, ftol=1e-3)
Compute the first TS step of the CS-TSSOS hierarchy for constrained complex polynomial optimization. If ConjugateBasis=true
, then include conjugate variables in monomial bases. If merge=true
, perform the PSD block merging. If solve=false
, then do not solve the SDP. If Gram=true
, then output the Gram matrix. If MomentOne=true
, add an extra first-order moment PSD constraint to the moment relaxation.
Input arguments
pop
: vector of the objective, inequality constraints, and equality constraintsz
: CPOP variables and their conjugated
: relaxation ordernb
: number of unit-norm variables inx
numeq
: number of equality constraintsCS
: method of chordal extension for correlative sparsity ("MF"
,"MD"
,false
)cliques
: the set of cliques used in correlative sparsityTS
: type of term sparsity ("block"
,"MD"
,"MF"
,false
)md
: tunable parameter for merging blocksnormality
: normal orderQUIET
: run in the quiet mode (true
,false
)rtol
: tolerance for rankgtol
: tolerance for global optimality gapftol
: tolerance for feasibility
Output arguments
opt
: optimumsol
: (near) optimal solution (ifsolution=true
)data
: other auxiliary data
opt,sol,data = complex_cs_tssos_first(supp::Vector{Vector{Vector{Vector{UInt16}}}}, coe::Vector{Vector{T}}, n, d;
nb=0, numeq=0, CS="MF", cliques=[], TS="block", merge=false, md=3, solver="Mosek", solution=false, dualize=false,
QUIET=false, solve=true, Gram=false, MomentOne=false, normality=!ConjugateBasis, cosmo_setting=cosmo_para(), mosek_setting=mosek_para(),
rtol=1e-2, gtol=1e-2, ftol=1e-3) where {T<:Number}
Compute the first TS step of the CS-TSSOS hierarchy for constrained complex polynomial optimization. Here the complex polynomial optimization problem is defined by supp
and coe
, corresponding to the supports and coeffients of pop
respectively.
TSSOS.complex_cs_tssos_higher!
— Functionopt,sol,data = complex_cs_tssos_higher!(data; TS="block", merge=false, md=3, QUIET=false, solve=true, dualize=false,
solution=false, Gram=false, MomentOne=false, cosmo_setting=cosmo_para(), mosek_setting=mosek_para())
Compute higher TS steps of the CS-TSSOS hierarchy.
TSSOS.add_complex_psatz!
— Functioninfo = add_complex_psatz!(model, nonneg, vars, ineq_cons, eq_cons, order; ipart=true, CS=false, cliques=[], TS="block",
SO=1, ConjugateBasis=false, normality=!ConjugateBasis, QUIET=false)
Add a complex Putinar's style Hermitian SOS representation of the complex polynomial nonneg
to the JuMP model
.
Input arguments
model
: a JuMP optimization modelnonneg
: a nonnegative complex polynomial constrained to be a complex Putinar's style Hermitian SOS on a semialgebraic setvars
: the set of POP variablesineq_cons
: inequality constraintseq_cons
: equality constraintsorder
: relaxation orderipart
: involving complex coefficientsCS
: method of chordal extension for correlative sparsity ("MF"
,"MD"
,"NC"
,false
)cliques
: the set of cliques used in correlative sparsityTS
: type of term sparsity ("block"
,"MD"
,"MF"
,false
)SO
: sparse orderConjugateBasis
: including conjugate variables in monomial bases.normality
: normal orderQUIET
: run in the quiet mode (true
,false
)
Output arguments
info
: auxiliary data