// Solve multiple right sides. If flags&UNIT then A diagonal is assumed to // to unit and is not referenced. (blas.TRSM) // alpha*B = A.-1*B if flags&LEFT // alpha*B = A.-T*B if flags&(LEFT|TRANS) // alpha*B = B*A.-1 if flags&RIGHT // alpha*B = B*A.-T if flags&(RIGHT|TRANS) // // Matrix A is N*N triangular matrix defined with flags bits as follow // LOWER non-unit lower triangular // LOWER|UNIT unit lower triangular // UPPER non-unit upper triangular // UPPER|UNIT unit upper triangular // // Matrix B is N*P if flags&LEFT or P*N if flags&RIGHT. // func SolveTrm(B, A *matrix.FloatMatrix, alpha float64, flags Flags) error { ok := true empty := false br, bc := B.Size() ar, ac := A.Size() switch flags & (LEFT | RIGHT) { case LEFT: empty = br == 0 ok = br == ac && ac == ar case RIGHT: empty = bc == 0 ok = bc == ar && ac == ar } if empty { return nil } if !ok { return onError("A, B size mismatch") } Ar := A.FloatArray() ldA := A.LeadingIndex() Br := B.FloatArray() ldB := B.LeadingIndex() E := bc if flags&RIGHT != 0 { E = br } // if more workers available can divide to tasks by B columns if flags&LEFT or by // B rows if flags&RIGHT. calgo.DSolveBlk(Br, Ar, alpha, calgo.Flags(flags), ldB, ldA, ac, 0, E, nB) return nil }
// Return Complex(Real, Imag). Return a new matrix. func Complex(Real, Imag *matrix.FloatMatrix) *matrix.ComplexMatrix { if !Real.SizeMatch(Imag.Size()) { return nil } C := matrix.ComplexZeros(Real.Size()) Rr := Real.FloatArray() Ir := Imag.FloatArray() Cr := C.ComplexArray() for i, _ := range Rr { Cr[i] = complex(Rr[i], Ir[i]) } return C }
// Generic matrix-matrix multpily. (blas.GEMM). Calculates // C = beta*C + alpha*A*B (default) // C = beta*C + alpha*A.T*B flags&TRANSA // C = beta*C + alpha*A*B.T flags&TRANSB // C = beta*C + alpha*A.T*B.T flags&(TRANSA|TRANSB) // // C is M*N, A is M*P or P*M if flags&TRANSA. B is P*N or N*P if flags&TRANSB. // func Mult(C, A, B *matrix.FloatMatrix, alpha, beta float64, flags Flags) error { var ok, empty bool // error checking must take in account flag values! ar, ac := A.Size() br, bc := B.Size() cr, cc := C.Size() switch flags & (TRANSA | TRANSB) { case TRANSA | TRANSB: empty = ac == 0 || br == 0 ok = cr == ac && cc == br && ar == bc case TRANSA: empty = ac == 0 || bc == 0 ok = cr == ac && cc == bc && ar == br case TRANSB: empty = ar == 0 || br == 0 ok = cr == ar && cc == br && ac == bc default: empty = ar == 0 || bc == 0 ok = cr == ar && cc == bc && ac == br } if empty { return nil } if !ok { return errors.New("Mult: size mismatch") } psize := int64(C.NumElements()) * int64(A.Cols()) Ar := A.FloatArray() ldA := A.LeadingIndex() Br := B.FloatArray() ldB := B.LeadingIndex() Cr := C.FloatArray() ldC := C.LeadingIndex() // matrix A, B common dimension P := A.Cols() if flags&TRANSA != 0 { P = A.Rows() } if nWorker <= 1 || psize <= limitOne { calgo.DMult(Cr, Ar, Br, alpha, beta, calgo.Flags(flags), ldC, ldA, ldB, P, 0, C.Cols(), 0, C.Rows(), vpLen, nB, mB) return nil } // here we have more than one worker available worker := func(cstart, cend, rstart, rend int, ready chan int) { calgo.DMult(Cr, Ar, Br, alpha, beta, calgo.Flags(flags), ldC, ldA, ldB, P, cstart, cend, rstart, rend, vpLen, nB, mB) ready <- 1 } colworks, rowworks := divideWork(C.Rows(), C.Cols(), nWorker) scheduleWork(colworks, rowworks, C.Cols(), C.Rows(), worker) return nil }
// Symmetric matrix multiply. (blas.SYMM) // C = beta*C + alpha*A*B (default) // C = beta*C + alpha*A.T*B flags&TRANSA // C = beta*C + alpha*A*B.T flags&TRANSB // C = beta*C + alpha*A.T*B.T flags&(TRANSA|TRANSB) // // C is N*P, A is N*N symmetric matrix. B is N*P or P*N if flags&TRANSB. // func MultSym(C, A, B *matrix.FloatMatrix, alpha, beta float64, flags Flags) error { var ok, empty bool ar, ac := A.Size() br, bc := B.Size() cr, cc := C.Size() switch flags & (TRANSA | TRANSB) { case TRANSA | TRANSB: empty = ac == 0 || br == 0 ok = ar == ac && cr == ac && cc == br && ar == bc case TRANSA: empty = ac == 0 || bc == 0 ok = ar == ac && cr == ac && cc == bc && ar == br case TRANSB: empty = ar == 0 || br == 0 ok = ar == ac && cr == ar && cc == br && ac == bc default: empty = ar == 0 || bc == 0 ok = ar == ac && cr == ar && cc == bc && ac == br } if empty { return nil } if !ok { return errors.New("MultSym: size mismatch") } /* if A.Rows() != A.Cols() { return errors.New("A matrix not square matrix."); } if A.Cols() != B.Rows() { return errors.New("A.cols != B.rows: size mismatch") } */ psize := int64(C.NumElements()) * int64(A.Cols()) Ar := A.FloatArray() ldA := A.LeadingIndex() Br := B.FloatArray() ldB := B.LeadingIndex() Cr := C.FloatArray() ldC := C.LeadingIndex() if nWorker <= 1 || psize <= limitOne { calgo.DMultSymm(Cr, Ar, Br, alpha, beta, calgo.Flags(flags), ldC, ldA, ldB, A.Cols(), 0, C.Cols(), 0, C.Rows(), vpLen, nB, mB) return nil } // here we have more than one worker available worker := func(cstart, cend, rstart, rend int, ready chan int) { calgo.DMultSymm(Cr, Ar, Br, alpha, beta, calgo.Flags(flags), ldC, ldA, ldB, A.Cols(), cstart, cend, rstart, rend, vpLen, nB, mB) ready <- 1 } colworks, rowworks := divideWork(C.Rows(), C.Cols(), nWorker) scheduleWork(colworks, rowworks, C.Cols(), C.Rows(), worker) return nil }
// Computes analytic center of A*x <= b with A m by n of rank n. // We assume that b > 0 and the feasible set is bounded. func Acent(A, b *matrix.FloatMatrix, niters int) (*matrix.FloatMatrix, []float64) { if niters <= 0 { niters = MAXITERS } ntdecrs := make([]float64, 0, niters) if A.Rows() != b.Rows() { return nil, nil } m, n := A.Size() x := matrix.FloatZeros(n, 1) H := matrix.FloatZeros(n, n) // Helper m*n matrix Dmn := matrix.FloatZeros(m, n) for i := 0; i < niters; i++ { // Gradient is g = A^T * (1.0/(b - A*x)). d = 1.0/(b - A*x) // d is m*1 matrix, g is n*1 matrix d := matrix.Minus(b, matrix.Times(A, x)).Inv() g := matrix.Times(A.Transpose(), d) // Hessian is H = A^T * diag(1./(b-A*x))^2 * A. // in the original python code expression d[:,n*[0]] creates // a m*n matrix where each column is copy of column 0. // We do it here manually. for i := 0; i < n; i++ { Dmn.SetColumn(i, d) } // Function mul creates element wise product of matrices. Asc := matrix.Mul(Dmn, A) blas.SyrkFloat(Asc, H, 1.0, 0.0, linalg.OptTrans) // Newton step is v = H^-1 * g. v := g.Copy().Scale(-1.0) lapack.PosvFloat(H, v) // Directional derivative and Newton decrement. lam := blas.DotFloat(g, v) ntdecrs = append(ntdecrs, math.Sqrt(-lam)) if ntdecrs[len(ntdecrs)-1] < TOL { fmt.Printf("last Newton decrement < TOL(%v)\n", TOL) return x, ntdecrs } // Backtracking line search. // y = d .* A*v y := d.Mul(A.Times(v)) step := 1.0 for 1-step*y.Max() < 0 { step *= BETA } search: for { // t = -step*y t := y.Copy().Scale(-step) // t = (1 + t) [e.g. t = 1 - step*y] t.Add(1.0) // ts = sum(log(1-step*y)) ts := t.Log().Sum() if -ts < ALPHA*step*lam { break search } step *= BETA } v.Scale(step) x = x.Plus(v) } // no solution !! fmt.Printf("Iteration %d exhausted\n", niters) return x, ntdecrs }
// Solves a convex optimization problem with a linear objective // // minimize c'*x // subject to f(x) <= 0 // G*x <= h // A*x = b. // // f is vector valued, convex and twice differentiable. The linear // inequalities are with respect to a cone C defined as the Cartesian // product of N + M + 1 cones: // // C = C_0 x C_1 x .... x C_N x C_{N+1} x ... x C_{N+M}. // // The first cone C_0 is the nonnegative orthant of dimension ml. The // next N cones are second order cones of dimension r[0], ..., r[N-1]. // The second order cone of dimension m is defined as // // { (u0, u1) in R x R^{m-1} | u0 >= ||u1||_2 }. // // The next M cones are positive semidefinite cones of order t[0], ..., t[M-1] >= 0. // // The structure of C is specified by DimensionSet dims which holds following sets // // dims.At("l") l, the dimension of the nonnegative orthant (array of length 1) // dims.At("q") r[0], ... r[N-1], list with the dimesions of the second-order cones // dims.At("s") t[0], ... t[M-1], array with the dimensions of the positive // semidefinite cones // // The default value for dims is l: []int{h.Rows()}, q: []int{}, s: []int{}. // // On exit Solution contains the result and information about the accurancy of the // solution. if SolutionStatus is Optimal then Solution.Result contains solutions // for the problems. // // Result.At("x")[0] primal solution // Result.At("snl")[0] non-linear constraint slacks // Result.At("sl")[0] linear constraint slacks // Result.At("y")[0] values for linear equality constraints y // Result.At("znl")[0] values of dual variables for nonlinear inequalities // Result.At("zl")[0] values of dual variables for linear inequalities // // If err is non-nil then sol is nil and err contains information about the argument or // computation error. // func Cpl(F ConvexProg, c, G, h, A, b *matrix.FloatMatrix, dims *sets.DimensionSet, solopts *SolverOptions) (sol *Solution, err error) { var mnl int var x0 *matrix.FloatMatrix mnl, x0, err = F.F0() if err != nil { return } if x0.Cols() != 1 { err = errors.New("'x0' must be matrix with one column") return } if c == nil { err = errors.New("'c' must be non nil matrix") return } if !c.SizeMatch(x0.Size()) { err = errors.New(fmt.Sprintf("'c' must be matrix of size (%d,1)", x0.Rows())) return } if h == nil { h = matrix.FloatZeros(0, 1) } if h.Cols() > 1 { err = errors.New("'h' must be matrix with 1 column") return } if dims == nil { dims = sets.NewDimensionSet("l", "q", "s") dims.Set("l", []int{h.Rows()}) } cdim := dims.Sum("l", "q") + dims.SumSquared("s") //cdim_pckd := dims.Sum("l", "q") + dims.SumPacked("s") //cdim_diag := dims.Sum("l", "q", "s") if h.Rows() != cdim { err = errors.New(fmt.Sprintf("'h' must be float matrix of size (%d,1)", cdim)) return } if G == nil { G = matrix.FloatZeros(0, c.Rows()) } if !G.SizeMatch(cdim, c.Rows()) { estr := fmt.Sprintf("'G' must be of size (%d,%d)", cdim, c.Rows()) err = errors.New(estr) return } // Check A and set defaults if it is nil if A == nil { // zeros rows reduces Gemv to vector products A = matrix.FloatZeros(0, c.Rows()) } if A.Cols() != c.Rows() { estr := fmt.Sprintf("'A' must have %d columns", c.Rows()) err = errors.New(estr) return } // Check b and set defaults if it is nil if b == nil { b = matrix.FloatZeros(0, 1) } if b.Cols() != 1 { estr := fmt.Sprintf("'b' must be a matrix with 1 column") err = errors.New(estr) return } if b.Rows() != A.Rows() { estr := fmt.Sprintf("'b' must have length %d", A.Rows()) err = errors.New(estr) return } var mc = matrixVar{c} var mb = matrixVar{b} var mA = matrixVarA{A} var mG = matrixVarG{G, dims} solvername := solopts.KKTSolverName if len(solvername) == 0 { if len(dims.At("q")) > 0 || len(dims.At("s")) > 0 { solvername = "chol" } else { solvername = "chol2" } } var factor kktFactor var kktsolver KKTCpSolver = nil if kktfunc, ok := solvers[solvername]; ok { // kkt function returns us problem spesific factor function. factor, err = kktfunc(G, dims, A, mnl) // solver is kktsolver = func(W *sets.FloatMatrixSet, x, z *matrix.FloatMatrix) (KKTFunc, error) { _, Df, H, err := F.F2(x, z) if err != nil { return nil, err } return factor(W, H, Df) } } else { err = errors.New(fmt.Sprintf("solver '%s' not known", solvername)) return } //return CplCustom(F, c, &mG, h, &mA, b, dims, kktsolver, solopts) return cpl_problem(F, &mc, &mG, h, &mA, &mb, dims, kktsolver, solopts, x0, mnl) }
// Solves a convex optimization problem with a linear objective // // minimize c'*x // subject to f(x) <= 0 // G*x <= h // A*x = b. // // using custom KTT equation solver and custom constraints G and A. // func CplCustomMatrix(F ConvexProg, c *matrix.FloatMatrix, G MatrixG, h *matrix.FloatMatrix, A MatrixA, b *matrix.FloatMatrix, dims *sets.DimensionSet, kktsolver KKTCpSolver, solopts *SolverOptions) (sol *Solution, err error) { var mnl int var x0 *matrix.FloatMatrix mnl, x0, err = F.F0() if err != nil { return } if x0.Cols() != 1 { err = errors.New("'x0' must be matrix with one column") return } if c == nil { err = errors.New("'c' must be non nil matrix") return } if !c.SizeMatch(x0.Size()) { err = errors.New(fmt.Sprintf("'c' must be matrix of size (%d,1)", x0.Rows())) return } if h == nil { h = matrix.FloatZeros(0, 1) } if h.Cols() > 1 { err = errors.New("'h' must be matrix with 1 column") return } if dims == nil { dims = sets.NewDimensionSet("l", "q", "s") dims.Set("l", []int{h.Rows()}) } cdim := dims.Sum("l", "q") + dims.SumSquared("s") if h.Rows() != cdim { err = errors.New(fmt.Sprintf("'h' must be float matrix of size (%d,1)", cdim)) return } // Check b and set defaults if it is nil if b == nil { b = matrix.FloatZeros(0, 1) } if b.Cols() != 1 { estr := fmt.Sprintf("'b' must be a matrix with 1 column") err = errors.New(estr) return } mc := matrixVar{c} mb := matrixVar{b} var mG MatrixVarG var mA MatrixVarA if G == nil { mG = &matrixVarG{matrix.FloatZeros(0, c.Rows()), dims} } else { mG = &matrixIfG{G} } if A == nil { mA = &matrixVarA{matrix.FloatZeros(0, c.Rows())} } else { mA = &matrixIfA{A} } return cpl_problem(F, &mc, mG, h, mA, &mb, dims, kktsolver, solopts, x0, mnl) }
// Solves a convex optimization problem with a linear objective // // minimize c'*x // subject to f(x) <= 0 // G*x <= h // A*x = b. // // using custom KTT equation solver. // func CplCustomKKT(F ConvexProg, c *matrix.FloatMatrix, G, h, A, b *matrix.FloatMatrix, dims *sets.DimensionSet, kktsolver KKTCpSolver, solopts *SolverOptions) (sol *Solution, err error) { var mnl int var x0 *matrix.FloatMatrix mnl, x0, err = F.F0() if err != nil { return } if x0.Cols() != 1 { err = errors.New("'x0' must be matrix with one column") return } if c == nil { err = errors.New("'c' must be non nil matrix") return } if !c.SizeMatch(x0.Size()) { err = errors.New(fmt.Sprintf("'c' must be matrix of size (%d,1)", x0.Rows())) return } if h == nil { h = matrix.FloatZeros(0, 1) } if h.Cols() > 1 { err = errors.New("'h' must be matrix with 1 column") return } if dims == nil { dims = sets.NewDimensionSet("l", "q", "s") dims.Set("l", []int{h.Rows()}) } cdim := dims.Sum("l", "q") + dims.SumSquared("s") if h.Rows() != cdim { err = errors.New(fmt.Sprintf("'h' must be float matrix of size (%d,1)", cdim)) return } if G == nil { G = matrix.FloatZeros(0, c.Rows()) } if !G.SizeMatch(cdim, c.Rows()) { estr := fmt.Sprintf("'G' must be of size (%d,%d)", cdim, c.Rows()) err = errors.New(estr) return } // Check A and set defaults if it is nil if A == nil { // zeros rows reduces Gemv to vector products A = matrix.FloatZeros(0, c.Rows()) } if A.Cols() != c.Rows() { estr := fmt.Sprintf("'A' must have %d columns", c.Rows()) err = errors.New(estr) return } // Check b and set defaults if it is nil if b == nil { b = matrix.FloatZeros(0, 1) } if b.Cols() != 1 { estr := fmt.Sprintf("'b' must be a matrix with 1 column") err = errors.New(estr) return } if b.Rows() != A.Rows() { estr := fmt.Sprintf("'b' must have length %d", A.Rows()) err = errors.New(estr) return } var mc = matrixVar{c} var mb = matrixVar{b} var mA = matrixVarA{A} var mG = matrixVarG{G, dims} return cpl_problem(F, &mc, &mG, h, &mA, &mb, dims, kktsolver, solopts, x0, mnl) }
func blockedQRT(A, T, W *matrix.FloatMatrix, nb int) { var ATL, ATR, ABL, ABR, AL, AR matrix.FloatMatrix var A00, A01, A02, A10, A11, A12, A20, A21, A22 matrix.FloatMatrix var WT, WB, W0, W1, W2 matrix.FloatMatrix var TTL, TTR, TBL, TBR matrix.FloatMatrix var T00, T01, T02, T11, T12, T22 matrix.FloatMatrix partition2x2( &ATL, &ATR, &ABL, &ABR, A, 0, 0, pTOPLEFT) partition2x2( &TTL, &TTR, &TBL, &TBR, T, 0, 0, pTOPLEFT) partition2x1( &WT, &WB, W, 0, pTOP) for ABR.Rows() > 0 && ABR.Cols() > 0 { repartition2x2to3x3(&ATL, &A00, &A01, &A02, &A10, &A11, &A12, &A20, &A21, &A22, A, nb, pBOTTOMRIGHT) repartition2x2to3x3(&TTL, &T00, &T01, &T02, nil, &T11, &T12, nil, nil, &T22, T, nb, pBOTTOMRIGHT) repartition2x1to3x1(&WT, &W0, &W1, &W2, W, nb, pBOTTOM) partition1x2( &AL, &AR, &ABR, nb, pLEFT) // current block size cb, rb := A11.Size() if rb < cb { cb = rb } // -------------------------------------------------------- // decompose left side AL == /A11\ // \A21/ unblockedQRT(&AL, &T11) // update A'tail i.e. A12 and A22 with (I - Y*T*Y.T).T * A'tail // compute: Q*T.C == C - Y*(C.T*Y*T).T updateWithQT(&A12, &A22, &A11, &A21, &T11, &W2, cb, true) // update T01: T01 = -T00*Y1.T*Y2*T11 // Y1 = /A10\ Y2 = /A11\ // \A20/ \A21/ // updateQRTReflector(&T01, &A10, &A20, &A11, &A21, &T00, &T11) // -------------------------------------------------------- continue3x3to2x2( &ATL, &ATR, &ABL, &ABR, &A00, &A11, &A22, A, pBOTTOMRIGHT) continue3x3to2x2( &TTL, &TTR, &TBL, &TBR, &T00, &T11, &T22, T, pBOTTOMRIGHT) continue3x1to2x1( &WT, &WB, &W0, &W1, W, pBOTTOM) } }
/* * Blocked QR decomposition with compact WY transform. As implemented * in lapack.xGEQRF subroutine. */ func blockedQR(A, Tvec, Twork, W *matrix.FloatMatrix, nb int) { var ATL, ATR, ABL, ABR, AL, AR matrix.FloatMatrix var A00, A01, A02, A10, A11, A12, A20, A21, A22 matrix.FloatMatrix var TT, TB matrix.FloatMatrix var t0, tau, t2, Tdiag, WT, WB, W0, W1, W2 matrix.FloatMatrix //var Twork, W *matrix.FloatMatrix Tdiag.DiagOf(Twork) partition2x2( &ATL, &ATR, &ABL, &ABR, A, 0, 0, pTOPLEFT) partition2x1( &TT, &TB, Tvec, 0, pTOP) partition2x1( &WT, &WB, W, 0, pTOP) for ABR.Rows() > 0 && ABR.Cols() > 0 { repartition2x2to3x3(&ATL, &A00, &A01, &A02, &A10, &A11, &A12, &A20, &A21, &A22, A, nb, pBOTTOMRIGHT) repartition2x1to3x1(&TT, &t0, &tau, &t2, Tvec, nb, pBOTTOM) repartition2x1to3x1(&WT, &W0, &W1, &W2, W, nb, pBOTTOM) partition1x2( &AL, &AR, &ABR, nb, pLEFT) // current block size cb, rb := A11.Size() if rb < cb { cb = rb } // -------------------------------------------------------- // decompose left side AL == /A11\ // \A21/ unblockedQRT(&AL, Twork) // copy scaling from T diagonal to tau-vector Tdiag.CopyTo(&tau) // update A'tail i.e. A12 and A22 with (I - Y*T*Y.T).T * A'tail // compute: C - Y*(C.T*Y*T).T updateWithQT(&A12, &A22, &A11, &A21, Twork, &W2, cb, true) // -------------------------------------------------------- continue3x3to2x2( &ATL, &ATR, &ABL, &ABR, &A00, &A11, &A22, A, pBOTTOMRIGHT) continue3x1to2x1( &TT, &TB, &t0, &tau, Tvec, pBOTTOM) continue3x1to2x1( &WT, &WB, &W0, &W1, W, pBOTTOM) } }
func qcl1(A, b *matrix.FloatMatrix) (*cvx.Solution, error) { // Returns the solution u, z of // // (primal) minimize || u ||_1 // subject to || A * u - b ||_2 <= 1 // // (dual) maximize b^T z - ||z||_2 // subject to || A'*z ||_inf <= 1. // // Exploits structure, assuming A is m by n with m >= n. m, n := A.Size() Fkkt := func(W *sets.FloatMatrixSet) (f cvx.KKTFunc, err error) { minor := 0 if !checkpnt.MinorEmpty() { minor = checkpnt.MinorTop() } err = nil f = nil beta := W.At("beta")[0].GetIndex(0) v := W.At("v")[0] // As = 2 * v *(v[1:].T * A) //v_1 := matrix.FloatNew(1, v.NumElements()-1, v.FloatArray()[1:]) v_1 := v.SubMatrix(1, 0).Transpose() As := matrix.Times(v, matrix.Times(v_1, A)).Scale(2.0) //As_1 := As.GetSubMatrix(1, 0, m, n) //As_1.Scale(-1.0) //As.SetSubMatrix(1, 0, matrix.Minus(As_1, A)) As_1 := As.SubMatrix(1, 0, m, n) As_1.Scale(-1.0) As_1.Minus(A) As.Scale(1.0 / beta) S := matrix.Times(As.Transpose(), As) checkpnt.AddMatrixVar("S", S) d1 := W.At("d")[0].SubMatrix(0, 0, n, 1).Copy() d2 := W.At("d")[0].SubMatrix(n, 0).Copy() // D = 4.0 * (d1**2 + d2**2)**-1 d := matrix.Plus(matrix.Mul(d1, d1), matrix.Mul(d2, d2)).Inv().Scale(4.0) // S[::n+1] += d S.Diag().Plus(d.Transpose()) err = lapack.Potrf(S) checkpnt.Check("00-Fkkt", minor) if err != nil { return } f = func(x, y, z *matrix.FloatMatrix) (err error) { minor := 0 if !checkpnt.MinorEmpty() { minor = checkpnt.MinorTop() } else { loopf += 1 minor = loopf } checkpnt.Check("00-f", minor) // -- z := - W**-T * z // z[:n] = -div( z[:n], d1 ) z_val := z.SubMatrix(0, 0, n, 1) z_res := matrix.Div(z_val, d1).Scale(-1.0) z.SubMatrix(0, 0, n, 1).Set(z_res) // z[n:2*n] = -div( z[n:2*n], d2 ) z_val = z.SubMatrix(n, 0, n, 1) z_res = matrix.Div(z_val, d2).Scale(-1.0) z.SubMatrix(n, 0, n, 1).Set(z_res) // z[2*n:] -= 2.0*v*( v[0]*z[2*n] - blas.dot(v[1:], z[2*n+1:]) ) v0_z2n := v.GetIndex(0) * z.GetIndex(2*n) v1_z2n := blas.DotFloat(v, z, &linalg.IOpt{"offsetx", 1}, &linalg.IOpt{"offsety", 2*n + 1}) z_res = matrix.Scale(v, -2.0*(v0_z2n-v1_z2n)) z.SubMatrix(2*n, 0, z_res.NumElements(), 1).Plus(z_res) // z[2*n+1:] *= -1.0 z.SubMatrix(2*n+1, 0).Scale(-1.0) // z[2*n:] /= beta z.SubMatrix(2*n, 0).Scale(1.0 / beta) // -- x := x - G' * W**-1 * z // z_n = z[:n], z_2n = z[n:2*n], z_m = z[-(m+1):], z_n := z.SubMatrix(0, 0, n, 1) z_2n := z.SubMatrix(n, 0, n, 1) z_m := z.SubMatrix(z.NumElements()-(m+1), 0) // x[:n] -= div(z[:n], d1) - div(z[n:2*n], d2) + As.T * z[-(m+1):] z_res = matrix.Minus(matrix.Div(z_n, d1), matrix.Div(z_2n, d2)) a_res := matrix.Times(As.Transpose(), z_m) z_res.Plus(a_res).Scale(-1.0) x.SubMatrix(0, 0, n, 1).Plus(z_res) // x[n:] += div(z[:n], d1) + div(z[n:2*n], d2) z_res = matrix.Plus(matrix.Div(z_n, d1), matrix.Div(z_2n, d2)) x.SubMatrix(n, 0, z_res.NumElements(), 1).Plus(z_res) checkpnt.Check("15-f", minor) // Solve for x[:n]: // // S*x[:n] = x[:n] - (W1**2 - W2**2)(W1**2 + W2**2)^-1 * x[n:] // w1 = (d1**2 - d2**2), w2 = (d1**2 + d2**2) w1 := matrix.Minus(matrix.Mul(d1, d1), matrix.Mul(d2, d2)) w2 := matrix.Plus(matrix.Mul(d1, d1), matrix.Mul(d2, d2)) // x[:n] += -mul( div(w1, w2), x[n:]) x_n := x.SubMatrix(n, 0) x_val := matrix.Mul(matrix.Div(w1, w2), x_n).Scale(-1.0) x.SubMatrix(0, 0, n, 1).Plus(x_val) checkpnt.Check("25-f", minor) // Solve for x[n:]: // // (d1**-2 + d2**-2) * x[n:] = x[n:] + (d1**-2 - d2**-2)*x[:n] err = lapack.Potrs(S, x) if err != nil { fmt.Printf("Potrs error: %s\n", err) } checkpnt.Check("30-f", minor) // Solve for x[n:]: // // (d1**-2 + d2**-2) * x[n:] = x[n:] + (d1**-2 - d2**-2)*x[:n] // w1 = (d1**-2 - d2**-2), w2 = (d1**-2 + d2**-2) w1 = matrix.Minus(matrix.Mul(d1, d1).Inv(), matrix.Mul(d2, d2).Inv()) w2 = matrix.Plus(matrix.Mul(d1, d1).Inv(), matrix.Mul(d2, d2).Inv()) x_n = x.SubMatrix(0, 0, n, 1) // x[n:] += mul( d1**-2 - d2**-2, x[:n]) x_val = matrix.Mul(w1, x_n) x.SubMatrix(n, 0, x_val.NumElements(), 1).Plus(x_val) checkpnt.Check("35-f", minor) // x[n:] = div( x[n:], d1**-2 + d2**-2) x_n = x.SubMatrix(n, 0) x_val = matrix.Div(x_n, w2) x.SubMatrix(n, 0, x_val.NumElements(), 1).Set(x_val) checkpnt.Check("40-f", minor) // x_n = x[:n], x-2n = x[n:2*n] x_n = x.SubMatrix(0, 0, n, 1) x_2n := x.SubMatrix(n, 0, n, 1) // z := z + W^-T * G*x // z[:n] += div( x[:n] - x[n:2*n], d1) x_val = matrix.Div(matrix.Minus(x_n, x_2n), d1) z.SubMatrix(0, 0, n, 1).Plus(x_val) checkpnt.Check("44-f", minor) // z[n:2*n] += div( -x[:n] - x[n:2*n], d2) x_val = matrix.Div(matrix.Plus(x_n, x_2n).Scale(-1.0), d2) z.SubMatrix(n, 0, n, 1).Plus(x_val) checkpnt.Check("48-f", minor) // z[2*n:] += As*x[:n] x_val = matrix.Times(As, x_n) z.SubMatrix(2*n, 0, x_val.NumElements(), 1).Plus(x_val) checkpnt.Check("50-f", minor) return nil } return } // matrix(n*[0.0] + n*[1.0]) c := matrix.FloatZeros(2*n, 1) c.SubMatrix(n, 0).SetIndexes(1.0) h := matrix.FloatZeros(2*n+m+1, 1) h.SetIndexes(1.0, 2*n) // h[2*n+1:] = -b h.SubMatrix(2*n+1, 0).Plus(b).Scale(-1.0) G := &matrixFs{A} dims := sets.DSetNew("l", "q", "s") dims.Set("l", []int{2 * n}) dims.Set("q", []int{m + 1}) var solopts cvx.SolverOptions solopts.ShowProgress = true if maxIter > 0 { solopts.MaxIter = maxIter } if len(solver) > 0 { solopts.KKTSolverName = solver } return cvx.ConeLpCustomMatrix(c, G, h, nil, nil, dims, Fkkt, &solopts, nil, nil) }