COSMO.ADMMIterationsType
ADMMIterations()

The number of ADMM iterations completed during the solve.

COSMO.AbstractComplexityWeightType
AbstractComplexityWeight

A supertype for weights that are based on the complexity of clique-related operations in the solver algorithm.

All subtypes have to define the function: compute_complexity_savings(n_1::Int, n_2::Int, n_m::Int, edge_weight::Subtype) --> Float64

with

  • n1: `|C1|`
  • n2: `|C2|`
  • n3: `|C1 ∪ C_2|`
COSMO.AbstractEdgeWeightType
AbstractEdgeWeight

A supertype for the weights that are put on the edges of the clique graph.

COSMO.AbstractMergeStrategyType
AbstractMergeStrategy

An abstract supertype for merge strategies.

A merge strategy determines how the cliques in a clique tree / clique graph are merged to improve computation time of the projection. Each merge strategy should implement the following functions:

  • initialise!: Initialise the graph / tree from the input clique graph / tree, allocate memory
  • traverse: A method that determines how the clique tree / graph is traversed
  • evaluate: A method that decides whether to merge two cliques
  • update_strategy!: A method to update local strategy-related information after a merge
COSMO.BoxType
Box(l, u)

Creates a box or intervall with lower boundary vector $l \in \mathbb{R}^m \cup \{-\infty\}^m$ and upper boundary vector$u \in \mathbb{R}^m\cup \{+\infty\}^m$.

COSMO.CliqueGraphMergeType
CliqueGraphMerge(edge_weight::AbstractEdgeWeight = ComplexityWeight()) <: AbstractGraphBasedMerge

The (default) merge strategy based on the reduced clique graph $\mathcal{G}(\mathcal{B}, \xi)$, for a set of cliques $\mathcal{B} = \{ \mathcal{C}_1, \dots, \mathcal{C}_p\}$ where the edge set $\xi$ is obtained by taking the edges of the union of clique trees.

Moreover, given an edge weighting function $e(\mathcal{C}_i,\mathcal{C}_j) = w_{ij}$, we compute a weight for each edge that quantifies the computational savings of merging the two cliques. After the initial weights are computed, we merge cliques in a loop:

while clique graph contains positive weights:

  • select two permissible cliques with the highest weight $w_{ij}$
  • merge cliques $\rightarrow$ update clique graph
  • recompute weights for updated clique graph

Custom edge weighting functions can be used by defining your own CustomEdgeWeight <: AbstractEdgeWeight and a corresponding edge_metric method. By default, the ComplexityWeight <: AbstractEdgeWeight is used which computes the weight based on the cardinalities of the cliques: $e(\mathcal{C}_i,\mathcal{C}_j) = |\mathcal{C}_i|^3 + |\mathcal{C}_j|^3 - |\mathcal{C}_i \cup \mathcal{C}_j|^3$.

See also: Garstka, Cannon, Goulart - A clique graph based merging strategy for decomposable SDPs (2019)

COSMO.ComplexityWeightType
ComplexityWeight

An edge weight that approximates the computational savings from merging two cliques Ci and Cj by e(Ci, Cj) = |Ci|^3 + |Cj|^3 - |Ci ∪ Cj|^3.

COSMO.ConstraintType
Constraint{T <: AbstractFloat}(A, b, convex_set_type, dim = 0, indices = 0:0)

Creates a COSMO constraint: Ax + b ∈ convex_set.

By default the following convex set types are supported: ZeroSet, Nonnegatives, SecondOrderCone, PsdCone, PsdConeTriangle.

Examples

julia> COSMO.Constraint([1 0;0 1], zeros(2), COSMO.Nonnegatives)
Constraint
Size of A: (2, 2)
ConvexSet: COSMO.Nonnegatives{Float64}

For convex sets that require their own data, it is possible to pass the pass the instantiated object directly rather than the type name.

Examples

julia> COSMO.Constraint([1 0;0 1], zeros(2), COSMO.Box([-1.;-1.],[1.;1.]))
Constraint
Size of A: (2, 2)
ConvexSet: COSMO.Box{Float64}

The optional arguments dim and indices can be used to specify A and b for subparts of variable x. If x has dimension dim = 4, then x[2] and x[3] can be constrained to the zero cone in the following way:

Examples

julia> c = COSMO.Constraint([1 0;0 1], zeros(2), COSMO.ZeroSet, 4, 2:3)
Constraint
Size of A: (2, 4)
ConvexSet: COSMO.ZeroSet{Float64}

Notice that extra columns of A have been added automatically.

julia>Matrix(c.A)
2×4 Array{Float64,2}:
0.0  1.0  0.0  0.0
0.0  0.0  1.0  0.0
COSMO.DualExponentialConeType
DualExponentialCone(MAX_ITERS::Int = 100, EXP_TOL = 1e-8)

Creates the dual exponential cone $\mathcal{K}^*_{exp} = \{(x, y, z) \mid x < 0, -xe^{y/x} \leq e^1 z \} \cup \{ (0,y,z) \mid y \geq 0, z \geq 0 \}$

COSMO.DualPowerConeType
DualPowerCone(alpha::Float64, MAX_ITERS::Int = 20, POW_TOL = 1e-8)

Creates the 3-d dual power cone $\mathcal{K}^*_{pow} = \{(u, v, w) \mid \left( \frac{u}{\alpha}\right)^\alpha \left( \frac{v}{1-\alpha}\right)^{(1-\alpha)} \geq \|w\|, u \geq 0, v \geq 0 \}$ with $0 < \alpha < 1$

COSMO.ExponentialConeType
ExponentialCone(MAX_ITERS = 100, EXP_TOL = 1e-8)

Creates the exponential cone $\mathcal{K}_{exp} = \{(x, y, z) \mid y \geq 0 ye^{x/y} ≤ z\} \cup \{ (x,y,z) \mid x \leq 0, y = 0, z \geq 0 \}$

COSMO.IterActivationType

Activate accelerator after start_iter iterations of the main algorithm.

COSMO.MergeLogType
	MergeLog

A struct to analyse the clique merges. Introduced for debugging purposes.

COSMO.ModelType
Model{T <: AbstractFloat}()

Initializes an empty COSMO model that can be filled with problem data using assemble!(model, P, q,constraints; [settings, x0, s0, y0]).

COSMO.NoMergeType
NoMerge <: AbstractMergeStrategy

A strategy that does not merge cliques.

COSMO.NonnegativesType
Nonnegatives(dim)

Creates the nonnegative orthant $\{ x \in \mathbb{R}^{dim} : x \ge 0 \}$ of dimension dim.

COSMO.ParentChildMergeType
ParentChildMerge(t_fill = 8, t_size = 8) <: AbstractTreeBasedMerge

The merge strategy suggested in Sun and Andersen - Decomposition in conic optimization with partially separable structure (2014). The initial clique tree is traversed in topological order and a clique $\mathcal{C}_\ell$ is greedily merged to its parent clique $\mathcal{C}_{par(\ell)}$ if at least one of the two conditions are met

  • $(| \mathcal{C}_{par(\ell)}| -| \eta_\ell|) (|\mathcal{C}_\ell| - |\eta_\ell|) \leq t_{\text{fill}}$ (fill-in condition)
  • $\max \left\{ |\nu_{\ell}|, |\nu_{par(\ell)}| \right\} \leq t_{\text{size}}$ (supernode size condition)
COSMO.PowerConeType
PowerCone(alpha::Float64, MAX_ITERS::Int = 20, POW_TOL = 1e-8)

Creates the 3-d power cone $\mathcal{K}_{pow} = \{(x, y, z) \mid x^\alpha y^{(1-\alpha)} \geq \|z\|, x \geq 0, y \geq 0 \}$ with $0 < \alpha < 1$

COSMO.PsdConeType
PsdCone(dim)

Creates the cone of symmetric positive semidefinite matrices $\mathcal{S}_+^{dim}$. The entries of the matrix X are stored column-by-column in the vector x of dimension dim. Accordingly $X \in \mathbb{S}_+ \Rightarrow x \in \mathcal{S}_+^{dim}$, where $X = \text{mat}(x)$.

COSMO.PsdConeTriangleType
PsdConeTriangle(dim)

Creates the cone of symmetric positive semidefinite matrices. The entries of the upper-triangular part of matrix X are stored in the vector x of dimension dim. A $r \times r$ matrix has $r(r+1)/2$ upper triangular elements and results in a vector of $\mathrm{dim} = r(r+1)/2$.

Examples

The matrix

\[\begin{bmatrix} x_1 & x_2 & x_4\\ x_2 & x_3 & x_5\\ x_4 & x_5 & x_6 \end{bmatrix}\]

is transformed to the vector $[x_1, x_2, x_3, x_4, x_5, x_6]^\top$ with corresponding constraint PsdConeTriangle(6).

COSMO.ResultType
Result{T <: AbstractFloat}

Object returned by the COSMO solver after calling optimize!(model). It has the following fields:

FieldnameTypeDescription
xVector{T}Primal variable
yVector{T}Dual variable
sVector{T}(Primal) set variable
obj_valTObjective value
iterIntTotal number of ADMM iterations (incl. safeguarding_iter)
safeguarding_iterIntNumber of iterations due to safeguarding of accelerator
statusSymbolSolution status
infoCOSMO.ResultInfoStruct with more information
timesCOSMO.ResultTimesStruct with several measured times
COSMO.ResultInfoType
ResultInfo{T <: AbstractFloat}

Object that contains further information about the primal residual, the dual residuals and the rho updates.

COSMO.ResultTimesType
ResultTimes

Part of the Result object returned by the solver. ResultTimes contains timing results for certain parts of the algorithm:

Time NameDescription
solver_timeTotal time used to solve the problem
setup_timeSetup time = graph_time + init_factor_time + scaling_time
scaling_timeTime to scale the problem data
graph_timeTime used to perform chordal decomposition
init_factor_timeTime used for initial factorisation of the system of linear equations
factor_update_timeSum of times used to refactor the system of linear equations due to rho
iter_timeTime spent in iteration loop
proj_timeTime spent in projection functions
post_timeTime used for post processing
update_timeTime spent in the update! function of the accelerator
accelerate_timeTime spent in the accelerate! function of the accelerator

By default COSMO only measures solver_time, setup_time and proj_time. To measure the other times set verbose_timing = true.

COSMO.SecondOrderConeType
SecondOrderCone(dim)

Creates the second-order cone (or Lorenz cone) $\{ (t,x) \in \mathrm{R}^{dim} : || x ||_2 \leq t \}$.

COSMO.SettingsType
COSMO.Settings{T}(; kwargs) where {T <: AbstractFloat}

Creates a COSMO settings object that is used to pass user settings to the solver.

ArgumentDescriptionValues (default)
rhoADMM rho step0.1
sigmaADMM sigma step1e-6
alphaRelaxation parameter1.6
eps_absAbsolute residual tolerance1e-5
eps_relRelative residual tolerance1e-5
nearly_ratioResidual tolerance ratio between MOI.NEARLY_FEASIBLE_POINT and $MOI.FEASIBLE_POINT$100
eps_prim_infPrimal infeasibility tolerance1e-5
eps_dual_infDual infeasibility tolerance1e-5
max_iterMaximum number of iterations5000
verboseVerbose printingfalse
verbose_timingVerbose timingfalse
kkt_solverLinear System solverQdldlKKTSolver
check_terminationCheck termination interval25
check_infeasibilityCheck infeasibility interval40
scalingNumber of scaling iterations10
adaptive_rhoAutomatic adaptation of step size parametertrue
adaptiverhomax_adaptionsMax number of rho adaptionstypemax(Int64) (deactivated)
decomposeActivate to decompose chordal psd constraintstrue
complete_dualActivate to complete the dual variable after decompositionfalse
merge_strategyChoose a strategy for clique mergingCliqueGraphMerge
compact_transformationChoose how a decomposed problem is transformedtrue
time_limitSet solver time limit in s0 (deactivated)
acceleratorAcceleration schemeAndersonAccelerator{T, Type2{QRDecomp}, RestartedMemory, NoRegularizer}
accelerator_activationAccelerator activationImmediateActivation
safeguardAccelerator safeguardingtrue
safeguard_tolSafeguarding tolerance2.0
COSMO.SuperNodeTreeType
SuperNodeTree

A structure to represent and analyse the sparsity pattern of the input matrix L.

Note:

Based on the merge_strategy in the constructor, SuperNodeTree might be initialised as a graph, i.e. seperators sep are left empty as well as snd_child and snd_post.

After cliques have been merged, a valid clique tree will be recomputed from the consolidated clique graph.

COSMO.WorkspaceType
Workspace{T <: AbstractFloat}()

Initializes an empty COSMO model that can be filled with problem data using assemble!(model, P, q,constraints; [settings, x0, s0, y0]).

COSMO.ZeroSetType
ZeroSet(dim)

Creates the zero set $\{ 0 \}^{dim}$ of dimension dim. If xZeroSet then all entries of x are zero.

COSMO._assemble_kkt_triangleMethod

Given P, A, sigma and rho return the upper / lower triangle, defined by shape, of the KKT condition matrix K.

COSMO.acceleration_post!Method
acceleration_post!(accelerator, ws, num_iter, safeguarding_iter)

A function that the accelerator can use after the nominal ADMM operator step.

In the case of an Anderson Accelerator this is used to check the quality of the accelerated candidate vector and take measures if the vector is of bad quality.

COSMO.acceleration_pre!Method
acceleration_pre!(accelerator, ws, num_iter)

A function that the accelerator can use before the nominal ADMM operator step.

In the case of an Anderson Accelerator this is used to calculate an accelerated candidate vector, that overwrites the current iterate.

COSMO.add_clique_entries!Method

Loop over all entries (i, j) in the clique and either set the correct row in A_I and b_I if (i, j) is not an overlap or add an overlap column with (-1 and +1) in the correct positions.

COSMO.admm_w!Method
admm_w!

ADMM-operator variable w update with over-relaxation parameter α.

COSMO.admm_x!Method
admm_x!

Evaluates the proximal operator of 1/2 x'Px + q'x s.t. Ax + s == b.

COSMO.admm_z!Method
admm_z!

Evaluates the proximal operator of the Indicator function I_{R^n × K}.

COSMO.allocate_sparse_matrixMethod

Given the row, column, and nzval vectors and dimensions, assemble the sparse matrix Aa of the decomposed problem in a slightly more memory efficient way.

COSMO.alternating_sequenceMethod

Returns the appropriate amount of memory for A.nzval, including, starting from n_start, the (+1 -1) entries for the overlaps.

COSMO.apply_rho_adaptation_rules!Method
apply_rho_adaptation_rules!()

Adapts the step-size parameter ρ when multiple conditions are satisfied.

  • If automatic rho interval is active, the adaptive_rho_interval is chosen as a fraction of the setup time.

  • Adapt rho if the iteration count reaches a multiple of adaptive_rho_interval and the maximum number of allowable updates has not been reached (adaptive_rho_max_adaptions).

  • If an accelerator is used, update rho at the next possible iterations, i.e. at the next non-accelerated iteration.

COSMO.assemble!Method
assemble!(model, P, q, constraint(s); [settings, x0, y0, s0])

Assembles a COSMO.Model with a cost function defind by P and q, and a number of constraints.

The positive semidefinite matrix P and vector q are used to specify the cost function of the optimization problem:

min   1/2 x'Px + q'x
s.t.  Ax + b ∈ C

constraints is a COSMO.Constraint or an array of COSMO.Constraint objects that are used to describe the constraints on x.


The optional keyword argument settings can be used to pass custom solver settings:

custom_settings = COSMO.Settings(verbose = true);
assemble!(model, P, q, constraints, settings = custom_settings)

The optional keyword arguments x0 and y0 can be used to provide the solver with warm starting values for the primal variable x and the dual variable y.

x_0 = [1.0; 5.0; 3.0]
COSMO.assemble!(model, P, q, constraints, x0 = x_0)
COSMO.augment_clique_based!Method
augment_clique_based!(ws::COSMO.Workspace{T})

If settings.compact_decomposition = true, decompose the problem based on the clique tree.

COSMO.check_termination!Method
check_termination!()

Checks the algorithms termination conditions at intervals specified in the settings and updates the algorithm status:

  • residuals are "small enough" –> :Solved
  • cost value out-ouf-range –> :Unsolved
  • infeasibility conditions satisfied –> :Primalinfeasible / :Dualinfeasible
  • time limit constraint reached –> :Timelimitreached
COSMO.clique_rows_mapMethod

Return the row ranges of each clique after the decomposition of C shifted by row_start.

COSMO.clique_tree_from_graph!Method
clique_tree_from_graph!(clique_graph::SuperNodeTree)

Given the cliques and edges of a clique graph, this function computes a valid clique tree.

This is necessary to perform the psd completion step of the dual variable after solving the problem.

COSMO.compute_reduced_clique_graph!Method
compute_reduced_clique_graph!(sep::Array{Set{Int}, 1}, snd::Array{Set{Int}, 1})

Compute the reduced clique graph (union of all clique trees) given an initial clique tree defined by its supernodes snd and separator sep sets.

We are using the algorithm described in Michel Habib and Juraj Stacho - Polynomial-time algorithm for the leafage ofchordal graphs (2009), which computes the reduced clique graph in the following way:

  1. Sort all minimal separators by size
  2. Initialise graph CG(R) with cliques as nodes and no edges
  3. for largest unprocessed separator S and | add an edge between any two cliques C1 and C2 if they both contain S and are in different connected components of CG(R) and store in edges. | Compute an edge weight used for merge decision and store in val. | Store the index of the separator which is the intersection C1 ∩ C2 in iter end
COSMO.compute_weights!Method

Compute the edge weight between all cliques specified by the edges (rows, cols).

COSMO.decomposed_dimMethod

For a cone C determine the sum of vector-dimensions of its cliques after decomposition (dim) and the total number of overlaps between them (overlaps).

COSMO.determine_parent_cliques!Method

Given the maximum weight spanning tree represented by E, determine a parent structure snd_par for the clique tree.

COSMO.edge_metricMethod
edge_metric(c_a::AbstractVector, c_b::AbstractVector, edge_weight::AbstractComplexityWeight)

Given two cliques c_a and c_b return a value for their edge weight.

COSMO.empty_model!Method
empty_model!(model)

Resets all the fields of model to that of a model created with COSMO.Model() (apart from the settings).

COSMO.evaluateMethod

Decide whether to merge the two cliques with clique indices cand.

COSMO.extra_columnsMethod

Returns the appropriate amount of memory for the columns of the augmented problem matrix A, including, starting from n_start, the columns for the (+1 -1) entries for the overlaps.

COSMO.fill_inMethod

Compute the amount of fill-in created by merging two cliques with the respective supernode and separator dimensions.

COSMO.find_A_dimensionMethod

Return the dimension of the problem after a clique tree based decomposition, given the sparsity patterns in sp_arr.

COSMO.find_componentsMethod

Find connected components in undirected separator graph represented by H.

COSMO.find_graph!Method
find_graph!(ci, rows::Array{Int, 1}, N::Int, C::AbstractConvexSet)

Given the indices of non-zero rows in rows:

  • Compute the sparsity pattern and find a chordal extension using QDLDL with AMD ordering F.perm.
  • If unconnected, connect the graph represented by the cholesky factor L
COSMO.find_neighborsMethod

find_neighbors(edges::SparseMatrixCSC, c::Int)

Find all the cliques connected to c which are given by the nonzeros in (c, 1:c-1) and (c+1:n, c).

COSMO.findnz!Method

Given a sparse matrix S, write the columns and non-zero values into the first numnz entries of J and V.

COSMO.get_blk_rowsMethod

For a PSD cone C of dimension nBlk, return the number of entries/rows in vectorized form.

COSMO.get_block_indicesMethod
get_block_indices(snd::Array{Int}, sep::Array{Int})

For a clique consisting of supernodes snd and seperators sep, compute all the indices (i, j) of the corresponding matrix block in the format (i, j, flag) where flag is equal to 0 if entry (i, j) corresponds to an overlap of the clique and 1 otherwise.

Nv is the number of vertices in the graph that we are trying to decompose.

COSMO.get_cliqueMethod

Return clique with post order ind (prevents returning empty arrays due to clique merging)

COSMO.get_decomposed_dimMethod

Returns the number of rows and the number of overlaps of all the blocks (cliques) represented in the tree after decomposition.

COSMO.get_nz_ind_mapMethod

A sparse vector that maps the index of svec(i, j) = k to the actual index of where that entry is stored in A.rowval.

COSMO.get_row_indexMethod

Given the svec index k and an offset row_range_col.start, return the location of the (i, j)th entry in the row vector rowval.

COSMO.isdisjointMethod
isdisjoint(c_a, c_b, is_sorted = false)

Checks whether two sets c_a, c_b have no common elements. Assumes by default that the sets are sorted.

COSMO.ispermissibleMethod

Check whether edge is permissible for a merge. An edge is permissible if for every common neighbor N, C1 ∩ N == C2 ∩ N or if no common neighbors exist.

COSMO.kruskal!Method
kruskal!(E::SparseMatrixCSC, num_cliques::Int)

Kruskal's algorithm to find a maximum weight spanning tree from the clique intersection graph.

E[i,j] holds the cardinalities of the intersection between two cliques (i, j). Changes the entries in the connectivity matrix E to a negative value if an edge between two cliques is included in the max spanning tree.

This is a modified version of https://github.com/JuliaGraphs/LightGraphs.jl/blob/master/src/spanningtrees/kruskal.jl

COSMO.log_merge!Method

Store information about the merge of the two merge candidates cand.

COSMO.mat_to_svec_indMethod

Convert the matrix indices (i, j) of A ∈ R^{m x n} into the corresponding linear index of v = svec(A, ::UpperTriangle).

COSMO.mat_to_vec_indMethod

Convert the matrix indices (i, j) of A ∈ R^{m x n} into the corresponding linear index of v = vec(A).

COSMO.max_elemMethod
max_elem(A::SparseMatrixCSC)

Find the matrix indices (i, j) of the first maximum element among the elements stored in A.nzval

COSMO.merge_child!Method

Merge two cliques cand of the tree t that are in a parent - child relationship.

COSMO.merge_cliques!Method

Given a supernode elimination tree t and a merge strategy strategy, merge the cliques in the tree.

COSMO.modify_clique_rows!Method

Given the nominal entry position k = svec(i, j) find and modify with new_row_val the actual location of that entry in the global row vector rowval.

COSMO.num_cone_decompositionMethod

Determine the total number of sets num_total after decomposition and the number of new psd cones num_new_psd_cones.

COSMO.optimize!Method

optimize!(model)

Attempts to solve the optimization problem defined in COSMO.Model object with the user settings defined in COSMO.Settings. Returns a COSMO.Result object.

COSMO.print_merge_logsMethod
print_merge_logs(ws; print_max = 10)

Print the merge logs for each decomposed cone in the problem.

COSMO.recover_μ!Method

The dual variable μ can be recovered from w, s via Moreau decomposition: μ = ρ (w - Π(w)).

COSMO.reorder_snd_consecutively!Method
	reorder_snd_consecutively!(sntree, ordering)

Takes a SuperNodeTree and reorders the vertices in each supernode (and separator) to have consecutive order.

The reordering is needed to achieve equal column structure for the psd completion of the dual variable Y. This also modifies ordering which maps the vertices in the sntree back to the actual location in the not reordered data, i.e. the primal constraint variable S and dual variables Y.

COSMO.reverse_decomposition!Method
reverse_decomposition!(ws::COSMO.Workspace, settings::COSMO.Settings)

After the problem has beend solved, undo the chordal decomposition.

Depending on the kind of transformation that was used, this involves:

  • Reassembling the original matrix S from its blocks
  • Reassembling the dual variable MU and performing a positive semidefinite completion.
COSMO.separator_graphMethod

Find the separator graph H given a separator and the relevant index-subset of cliques.

COSMO.set!Method
set!(model, P, q, A, b, convex_sets, [settings])

Sets model data directly based on provided fields.

COSMO.sort_col_wise!Method

Given the row vector and colptr of a sparse matrix, sort the entries within rowvals and nzvals withing the first n0 columns.

COSMO.split_cliques!Method

Traverse the clique tree in descending topological order and split the clique sets into supernodes and separators.

COSMO.symmetrize_full!Method
symmetrize_full!(A)

Symmetrizes the matrix A by calculating A = 0.5 * (A + A') and storing the result in-place.

COSMO.symmetrize_upper!Method
symmetrize_upper!(A)

Symmetrizes the matrix A by calculating A = 0.5 * (A + A') but only performs the operation on the upper triangular part.

COSMO.traverseMethod

Find the next two cliques in the clique graph t to merge.

COSMO.traverseMethod

Traverse tree t in descending topological order and return parent and clique (root has highest order).

COSMO.type_checksMethod

Check whether the model will contain any PSD constraints with unsupported Floating-point precision.

COSMO.union_dimMethod

Find the size of the set A ∪ B under the assumption that A and B only have unique elements.

COSMO.update!Method
update!(model, q = nothing, b = nothing)

Updates the model data for b or q. This can be done without refactoring the KKT matrix. The vectors will be appropriatly scaled.

COSMO.vec_dimMethod

Given the side dimension of a PSD cone return the number of stored entries.

COSMO.warm_start!Method
warm_start!(model, x0, y0)

Provides the COSMO.Model with warm starting values for the primal variable x and the dual variable y.

COSMO.warm_start_dual!Method
warm_start_dual!(model, y0, [ind])

Provides the COSMO.Model with warm starting values for the dual variable y. ind can be used to warm start certain components of y.

COSMO.warm_start_primal!Method
warm_start_primal!(model, x0, [ind])

Provides the COSMO.Model with warm starting values for the primal variable x. ind can be used to warm start certain components of x.

COSMO.warm_start_slack!Method
warm_start_slack!(model, s0, [ind])

Provides the COSMO.Model with warm starting values for the primal slack variable s. ind can be used to warm start certain components of s.