utils.func_from_others.Fawzi_func¶
- utils.func_from_others.Fawzi_func.matrix_geo_mean_hypo_cone(sz, t, iscplx, fullhyp)¶
- MATRIX_GEO_MEAN_HYPO_CONE Matrix geometric mean cone.
MATRIX_GEO_MEAN_HYPO_CONE(sz,t,iscplx,1) returns a CVX tuple {A,B,T} of matrices constrained to satisfy
A #_{t} B >= T
- where:
A #_{t} B is the t-weighted geometric mean of A and B (see
reference below for details). Parameter t should be in [0,1]. * The inequality “>=” is in the positive semidefinite order * A,B,T are symmetric (Hermitian if iscplx=1) matrices of size sz
- Note on parameter sz:
It is possible to provide for sz an array [sz(1) sz(2) sz(3) …]. This will return an array of prod(sz(3:end)) triples, where each triple has size sz(1)=sz(2) and is constrained to live in the operator relative entropy ocne. Note that sz(1) and sz(2) must be equal in this case. See documentation for CVX’s “sets/semidefinite.m” for more information.
- Note on parameter fullhyp:
- In many applications one doesn’t need the full hypograph
hyp_t = {(A,B,T) : A #_{t} B >= T}
- but rather it is enough to work with a convex set C_t that satisfies
(A,B,A #_{t} B) in C_t (A,B,T) in C_t => A #_{t} B >= T
In this case one should set fullhyp = 0. The SDP description will be (slightly) smaller. (By default fullhyp is set to 1).
- AUTHORS
Hamza Fawzi and James Saunderson
- REFERENCE
This code is based on the paper: “Lieb’s concavity theorem, matrix geometric means and semidefinite optimization” by Hamza Fawzi and James Saunderson (arXiv:1512.03401)
- utils.func_from_others.Fawzi_func.op_rel_entr_epi_cone(sz, iscplx, m, k, e, apx)¶
- OP_REL_ENTR_EPI_CONE Operator relative entropy cone
OP_REL_ENTR_EPI_CONE(sz,iscplx) returns a CVX triple {X,Y,TAU} of matrices constrained to satisfy
X^{1/2} * logm(X^{1/2}*Y^{-1}*X^{1/2}) * X^{1/2} <= TAU
- where:
The inequality “<=” is in the positive semidefinite order
X,Y,TAU are symmetric (Hermitian if iscplx=1) matrices of size sz
logm is the matrix logarithm (in base e)
- Usage example:
variable X(4,4) symmetric variable Y(4,4) symmetric variable TAU(4,4) symmetric {X,Y,TAU} == op_rel_entr_epi_cone(4,0)
- Note on parameter sz:
It is possible to provide for sz an array [sz(1) sz(2) sz(3) …]. This will return an array of prod(sz(3:end)) triples, where each triple has size sz(1)=sz(2) and is constrained to live in the operator relative entropy ocne. Note that sz(1) and sz(2) must be equal in this case. See documentation for CVX’s “sets/semidefinite.m” for more information.
OP_REL_ENTR_EPI_CONE(sz,iscplx,m,k) The additional parameters m and k can be used to control the accuracy of the approximation (see reference for more details). m is the number of quadrature nodes to use and k is the number of square-roots to take. Default (m,k) = (3,3).
OP_REL_ENTR_EPI_CONE(sz,iscplx,m,k,e) If an additional matrix e is provided of size nxr (where n=sz(1)), the returned tuple(s) {X,Y,TAU} is constrained to satisfy instead:
e’*(D_{op}(X||Y))*e <= TAU.
Note that TAU here is of size rxr. The default case corresponds to e=eye(n). When r is small compared to n this can be helpful to reduce the size of small LMIs from 2nx2n to (n+r)x(n+r).
OP_REL_ENTR_EPI_CONE(sz,iscplx,m,k,e,apx) The last parameter apx indicates which approximation of the logarithm to use. Possible values: - apx = +1: Upper approximation on D_{op} [inner approximation of the
operator relative entropy cone]
- apx = -1: Lower approximation on D_{op} [outer approximation of the
operator relative entropy cone]
- apx = 0 (Default): Pade approximation (neither upper nor lower) but
slightly better accuracy than apx = +1 or -1.
- AUTHORS
Hamza Fawzi, James Saunderson and Pablo A. Parrilo
- REFERENCE
This code is based on the paper: “Semidefinite approximations of matrix logarithm” by Hamza Fawzi, James Saunderson and Pablo A. Parrilo
- utils.func_from_others.Fawzi_func.quantum_entr(X, m, k, apx)¶
- QUANTUM_ENTR Quantum (von Neumann) entropy
QUANTUM_ENTR(X) returns -trace(X*logm(X)) where the logarithm is base e (and not base 2!). X must be a positive semidefinite matrix.
- Disciplined convex programming information:
QUANTUM_ENTR is concave in X. This function implements the semidefinite programming approximation given in the reference below.
Parameters m and k control the accuracy of this approximation: m is the number of quadrature nodes to use and k the number of square-roots to take. See reference for more details. Default (m,k) = (3,3).
Parameter apx indicates which approximation r of the von Neumann entropy to use: - apx = +1: Upper approximation of entropy (H(X) <= r(X)) - apx = -1: Lower approximation of entropy (r(X) <= H(X)) - apx = 0 (Default): Pade approximation (neither upper nor lower),
but slightly better accuracy than apx=+1 or -1.
The upper and lower approximation are based on rational functions derived from Gauss-Radau quadrature, see documentation in the ‘doc’ folder.
REQUIRES: op_rel_entr_epi_cone Implementation uses H(X) = -trace( D_{op}(X||I) ) where D_{op} is the operator relative entropy:
D_{op}(X||Y) = X^{1/2}*logm(X^{1/2} Y^{-1} X^{1/2})*X^{1/2}
- AUTHORS
Hamza Fawzi, James Saunderson and Pablo A. Parrilo
- REFERENCE
This code is based on the paper: “Semidefinite approximations of matrix logarithm” by Hamza Fawzi, James Saunderson and Pablo A. Parrilo
- utils.func_from_others.Fawzi_func.quantum_rel_entr(A, B, m, k, apx)¶
- QUANTUM_REL_ENTR Quantum relative entropy
QUANTUM_REL_ENTR(A,B) returns trace(A*(logm(A)-logm(B))) where A and B are positive semidefinite matrices such that im(A) subseteq im(B) (otherwise the function evaluates to infinity). Note this function uses logarithm base e (and not base 2!).
- Disciplined convex programming information:
QUANTUM_REL_ENTR(A,B) is convex (jointly) in (A,B) This function implements the semidefinite programming approximation given in the reference below.
Parameters m and k control the accuracy of this approximation: m is the number of quadrature nodes to use and k the number of square-roots to take. See reference for more details. Default (m,k) = (3,3).
Parameter apx indicates which approximation r of the relative entropy function to use: - apx = +1: Upper approximation (D(A|B) <= r(A,B)) - apx = -1: Lower approximation (r(A,B) <= D(A|B)) - apx = 0 (Default): Pade approximation (neither upper nor lower),
but slightly better accuracy than apx=+1 or -1.
The upper and lower approximation are based on rational functions derived from Gauss-Radau quadrature, see documentation in the ‘doc’ folder.
REQUIRES: op_rel_entr_epi_cone Implementation uses the expression
D(A||B) = e’*D_{op} (A otimes I || I otimes B) )*e
where D_{op} is the operator relative entropy and e = In(:) where In is the nxn identity matrix.
- AUTHORS
Hamza Fawzi, James Saunderson and Pablo A. Parrilo
- REFERENCE
This code is based on the paper: “Semidefinite approximations of matrix logarithm” by Hamza Fawzi, James Saunderson and Pablo A. Parrilo
- utils.func_from_others.Fawzi_func.trace_logm(X, C, m, k, apx)¶
- TRACE_LOGM Trace of logarithm
TRACE_LOGM(X) returns trace(logm(X)) where X is a positive definite matrix. TRACE_LOGM(X,C) returns trace(C*logm(X)) where C is a positive semidefinite matrix of the same size as X.
- Disciplined convex programming information:
TRACE_LOGM(X,C) is concave in X, provided C is a (fixed) positive semidefinite matrix. This function implements the semidefinite programming approximation given in the reference below.
Parameters m and k control the accuracy of this approximation: m is the number of quadrature nodes to use and k the number of square-roots to take. See reference for more details. Default (m,k) = (3,3).
Parameter apx indicates which approximation r of logm(X) to use: - apx = +1: Upper approximation (logm(X) <= r(X)) - apx = -1: Lower approximation (r(X) <= logm(X)) - apx = 0 (Default): Pade approximation (neither upper nor lower),
but slightly better accuracy than apx=+1 or -1.
The upper and lower approximation are based on rational functions derived from Gauss-Radau quadrature, see documentation in the ‘doc’ folder.
REQUIRES: op_rel_entr_epi_cone Implementation uses trace(C*logm(X)) = -trace(C*D_{op}(I||X)) where D_{op} is the operator relative entropy:
D_{op}(X||Y) = X^{1/2}*logm(X^{1/2} Y^{-1} X^{1/2})*X^{1/2}
- AUTHORS
Hamza Fawzi, James Saunderson and Pablo A. Parrilo
- REFERENCE
This code is based on the paper: “Semidefinite approximations of matrix logarithm” by Hamza Fawzi, James Saunderson and Pablo A. Parrilo