func TestBidiagReduceUnblocked(t *testing.T) { N := 217 M := 269 conf := gomas.NewConf() conf.LB = 0 zeromean := cmat.NewFloatNormSource() A := cmat.NewMatrix(M, N) A.SetFrom(zeromean) tauq := cmat.NewMatrix(N, 1) taup := cmat.NewMatrix(N, 1) At := cmat.NewMatrix(N, M) blasd.Transpose(At, A) tauqt := cmat.NewMatrix(N, 1) taupt := cmat.NewMatrix(N, 1) W := lapackd.Workspace(M + N) lapackd.BDReduce(A, tauq, taup, W, conf) lapackd.BDReduce(At, tauqt, taupt, W, conf) // BiRed(A) == BiRed(A.T).T blasd.Plus(At, A, 1.0, -1.0, gomas.TRANSB) blasd.Axpy(tauqt, taup, -1.0) blasd.Axpy(taupt, tauq, -1.0) nrm := lapackd.NormP(At, lapackd.NORM_ONE) t.Logf("M=%d, N=%d || BiRed(A) - BiRed(A.T).T||_1 : %e\n", M, N, nrm) nrm = lapackd.NormP(taupt, lapackd.NORM_ONE) t.Logf(" || BiRed(A).tauq - BiRed(A.T).taup||_1 : %e\n", nrm) nrm = lapackd.NormP(tauqt, lapackd.NORM_ONE) t.Logf(" || BiRed(A).taup - BiRed(A.T).tauq||_1 : %e\n", nrm) }
func TestTriRedUpper(t *testing.T) { N := 843 nb := 48 conf := gomas.NewConf() conf.LB = 0 A := cmat.NewMatrix(N, N) tau := cmat.NewMatrix(N, 1) src := cmat.NewFloatNormSource() A.SetFrom(src, cmat.UPPER) A1 := cmat.NewCopy(A) tau1 := cmat.NewCopy(tau) _ = A1 W := lapackd.Workspace(N) W1 := lapackd.Workspace(N * nb) lapackd.TRDReduce(A, tau, W, gomas.UPPER, conf) conf.LB = nb lapackd.TRDReduce(A1, tau1, W1, gomas.UPPER, conf) blasd.Plus(A, A1, -1.0, 1.0, gomas.NONE) nrm := lapackd.NormP(A, lapackd.NORM_ONE) t.Logf("N=%d, ||unblk.Trired(A) - blk.Trired(A)||_1: %e\n", N, nrm) blasd.Axpy(tau, tau1, -1.0) nrm = blasd.Nrm2(tau) t.Logf(" ||unblk.Trired(tau) - blk.Trired(tau)||_1: %e\n", nrm) }
/* * Applies a real elementary reflector H to a real m by n matrix A, * from either the left or the right. H is represented in the form * * H = I - tau * ( 1 ) * ( 1 v.T ) * ( v ) * * where tau is a real scalar and v is a real vector. * * If tau = 0, then H is taken to be the unit cmat. * * A is /a1\ a1 := a1 - w1 * \A2/ A2 := A2 - v*w1 * w1 := tau*(a1 + A2.T*v) if side == LEFT * := tau*(a1 + A2*v) if side == RIGHT * * Intermediate work space w1 required as parameter, no allocation. */ func applyHouseholder2x1(tau, v, a1, A2, w1 *cmat.FloatMatrix, flags int) *gomas.Error { var err *gomas.Error = nil tval := tau.Get(0, 0) if tval == 0.0 { return err } // shape oblivious vector copy. blasd.Axpby(w1, a1, 1.0, 0.0) if flags&gomas.LEFT != 0 { // w1 = a1 + A2.T*v err = blasd.MVMult(w1, A2, v, 1.0, 1.0, gomas.TRANSA) } else { // w1 = a1 + A2*v err = blasd.MVMult(w1, A2, v, 1.0, 1.0, gomas.NONE) } // w1 = tau*w1 blasd.Scale(w1, tval) // a1 = a1 - w1 blasd.Axpy(a1, w1, -1.0) // A2 = A2 - v*w1 if flags&gomas.LEFT != 0 { err = blasd.MVUpdate(A2, v, w1, -1.0) } else { err = blasd.MVUpdate(A2, w1, v, -1.0) } return err }
/* * Reduce upper triangular matrix to tridiagonal. * * Elementary reflectors Q = H(n-1)...H(2)H(1) are stored on upper * triangular part of A. Reflector H(n-1) saved at column A(n) and * scalar multiplier to tau[n-1]. If parameter `tail` is true then * this function is used to reduce tail part of partially reduced * matrix and tau-vector partitioning is starting from last position. */ func unblkReduceTridiagUpper(A, tauq, W *cmat.FloatMatrix, tail bool) { var ATL, ABR cmat.FloatMatrix var A00, a01, a11, A22 cmat.FloatMatrix var tqT, tqB, tq0, tauq1, tq2 cmat.FloatMatrix var y21 cmat.FloatMatrix var v0 float64 toff := 1 if tail { toff = 0 } util.Partition2x2( &ATL, nil, nil, &ABR, A, 0, 0, util.PBOTTOMRIGHT) util.Partition2x1( &tqT, &tqB, tauq, toff, util.PBOTTOM) for n(&ATL) > 0 { util.Repartition2x2to3x3(&ATL, &A00, &a01, nil, nil, &a11, nil, nil, nil, &A22, A, 1, util.PTOPLEFT) util.Repartition2x1to3x1(&tqT, &tq0, &tauq1, &tq2, tauq, 1, util.PTOP) // set temp vectors for this round y21.SetBuf(n(&A00), 1, n(&A00), W.Data()) // ------------------------------------------------------ // Compute householder to zero super-diagonal entries computeHouseholderRev(&a01, &tauq1) tauqv := tauq1.Get(0, 0) // set superdiagonal to unit v0 = a01.Get(-1, 0) a01.Set(-1, 0, 1.0) // y21 := A22*a12t blasd.MVMultSym(&y21, &A00, &a01, tauqv, 0.0, gomas.UPPER) // beta := tauq*a12t*y21 beta := tauqv * blasd.Dot(&a01, &y21) // y21 := y21 - 0.5*beta*a125 blasd.Axpy(&y21, &a01, -0.5*beta) // A22 := A22 - a12t*y21.T - y21*a12.T blasd.MVUpdate2Sym(&A00, &a01, &y21, -1.0, gomas.UPPER) // restore superdiagonal value a01.Set(-1, 0, v0) // ------------------------------------------------------ util.Continue3x3to2x2( &ATL, nil, nil, &ABR, &A00, &a11, &A22, A, util.PTOPLEFT) util.Continue3x1to2x1( &tqT, &tqB, &tq0, &tauq1, tauq, util.PTOP) } }
func test_bdsvd(N, flags, kind int, verbose bool, t *testing.T) { var At, sD, sE, tmp cmat.FloatMatrix uplo := "upper" offdiag := 1 if flags&gomas.LOWER != 0 { offdiag = -1 uplo = "lower" } A0 := cmat.NewMatrix(N, N) desc := setDiagonals(A0, offdiag, kind) At.SubMatrix(A0, 0, 0, N, N) sD.Diag(A0, 0) sE.Diag(A0, offdiag) D := cmat.NewCopy(&sD) E := cmat.NewCopy(&sE) // unit singular vectors U := cmat.NewMatrix(N, N) sD.Diag(U, 0) sD.Add(1.0) V := cmat.NewMatrix(N, N) sD.Diag(V, 0) sD.Add(1.0) W := cmat.NewMatrix(4*N, 1) C := cmat.NewMatrix(N, N) lapackd.BDSvd(D, E, U, V, W, flags|gomas.WANTU|gomas.WANTV) blasd.Mult(C, U, U, 1.0, 0.0, gomas.TRANSA) sD.Diag(C) sD.Add(-1.0) nrmu := lapackd.NormP(C, lapackd.NORM_ONE) blasd.Mult(C, V, V, 1.0, 0.0, gomas.TRANSB) sD.Add(-1.0) nrmv := lapackd.NormP(C, lapackd.NORM_ONE) blasd.Mult(C, U, A0, 1.0, 0.0, gomas.TRANSA) blasd.Mult(&At, C, V, 1.0, 0.0, gomas.TRANSB) if verbose && N < 10 { t.Logf("D:\n%v\n", asRow(&tmp, D)) t.Logf("U:\n%v\n", U) t.Logf("V:\n%v\n", V) t.Logf("U.T*A*V\n%v\n", &At) } sD.Diag(&At) blasd.Axpy(&sD, D, -1.0) nrma := lapackd.NormP(&At, lapackd.NORM_ONE) t.Logf("N=%d [%s,%s] ||U.T*A*V - bdsvd(A)||_1: %e\n", N, uplo, desc, nrma) t.Logf(" ||I - U.T*U||_1: %e\n", nrmu) t.Logf(" ||I - V.T*V||_1: %e\n", nrmv) }
/* * Tridiagonal reduction of LOWER triangular symmetric matrix, zero elements below 1st * subdiagonal: * * A = (1 - tau*u*u.t)*A*(1 - tau*u*u.T) * = (I - tau*( 0 0 )) (a11 a12) (I - tau*( 0 0 )) * ( ( 0 u*u.t)) (a21 A22) ( ( 0 u*u.t)) * * a11, a12, a21 not affected * * from LEFT: * A22 = A22 - tau*u*u.T*A22 * from RIGHT: * A22 = A22 - tau*A22*u.u.T * * LEFT and RIGHT: * A22 = A22 - tau*u*u.T*A22 - tau*(A22 - tau*u*u.T*A22)*u*u.T * = A22 - tau*u*u.T*A22 - tau*A22*u*u.T + tau*tau*u*u.T*A22*u*u.T * [x = tau*A22*u (vector)] (SYMV) * A22 = A22 - u*x.T - x*u.T + tau*u*u.T*x*u.T * [beta = tau*u.T*x (scalar)] (DOT) * = A22 - u*x.T - x*u.T + beta*u*u.T * = A22 - u*(x - 0.5*beta*u).T - (x - 0.5*beta*u)*u.T * [w = x - 0.5*beta*u] (AXPY) * = A22 - u*w.T - w*u.T (SYR2) * * Result of reduction for N = 5: * ( d . . . . ) * ( e d . . . ) * ( v1 e d . . ) * ( v1 v2 e d . ) * ( v1 v2 v3 e d ) */ func unblkReduceTridiagLower(A, tauq, W *cmat.FloatMatrix) { var ATL, ABR cmat.FloatMatrix var A00, a11, a21, A22 cmat.FloatMatrix var tqT, tqB, tq0, tauq1, tq2 cmat.FloatMatrix var y21 cmat.FloatMatrix var v0 float64 util.Partition2x2( &ATL, nil, nil, &ABR, A, 0, 0, util.PTOPLEFT) util.Partition2x1( &tqT, &tqB, tauq, 0, util.PTOP) for m(&ABR) > 0 && n(&ABR) > 0 { util.Repartition2x2to3x3(&ATL, &A00, nil, nil, nil, &a11, nil, nil, &a21, &A22, A, 1, util.PBOTTOMRIGHT) util.Repartition2x1to3x1(&tqT, &tq0, &tauq1, &tq2, tauq, 1, util.PBOTTOM) // set temp vectors for this round y21.SetBuf(n(&A22), 1, n(&A22), W.Data()) // ------------------------------------------------------ // Compute householder to zero subdiagonal entries computeHouseholderVec(&a21, &tauq1) tauqv := tauq1.Get(0, 0) // set subdiagonal to unit v0 = a21.Get(0, 0) a21.Set(0, 0, 1.0) // y21 := tauq*A22*a21 blasd.MVMultSym(&y21, &A22, &a21, tauqv, 0.0, gomas.LOWER) // beta := tauq*a21.T*y21 beta := tauqv * blasd.Dot(&a21, &y21) // y21 := y21 - 0.5*beta*a21 blasd.Axpy(&y21, &a21, -0.5*beta) // A22 := A22 - a21*y21.T - y21*a21.T blasd.MVUpdate2Sym(&A22, &a21, &y21, -1.0, gomas.LOWER) // restore subdiagonal a21.Set(0, 0, v0) // ------------------------------------------------------ util.Continue3x3to2x2( &ATL, nil, nil, &ABR, &A00, &a11, &A22, A, util.PBOTTOMRIGHT) util.Continue3x1to2x1( &tqT, &tqB, &tq0, &tauq1, tauq, util.PBOTTOM) } }
func TestReduceBidiagBlkWide(t *testing.T) { N := 911 M := 823 nb := 48 conf := gomas.NewConf() conf.LB = 0 zeromean := cmat.NewFloatNormSource() A := cmat.NewMatrix(M, N) A.SetFrom(zeromean) tauq := cmat.NewMatrix(N, 1) taup := cmat.NewMatrix(N, 1) A1 := cmat.NewCopy(A) tauq1 := cmat.NewMatrix(N, 1) taup1 := cmat.NewMatrix(N, 1) W := lapackd.Workspace(M + N) lapackd.BDReduce(A, tauq, taup, W, conf) conf.LB = nb W1 := lapackd.Workspace(lapackd.BDReduceWork(A1, conf)) lapackd.BDReduce(A1, tauq1, taup1, W1, conf) // BiRed(A) == BiRed(A.T).T blasd.Plus(A1, A, 1.0, -1.0, gomas.NONE) blasd.Axpy(tauq1, tauq, -1.0) blasd.Axpy(taup1, taup, -1.0) nrm := lapackd.NormP(A1, lapackd.NORM_ONE) t.Logf("M=%d, N=%d || BiRed(A) - blk.BiRed(A)||_1 : %e\n", M, N, nrm) nrm = lapackd.NormP(taup1, lapackd.NORM_ONE) t.Logf(" || BiRed(A).tauq - blk.BiRed(A).taup||_1 : %e\n", nrm) nrm = lapackd.NormP(tauq1, lapackd.NORM_ONE) t.Logf(" || BiRed(A).taup - blk.BiRed(A).tauq||_1 : %e\n", nrm) }
func TestBidiagReduceBlockedTall(t *testing.T) { N := 711 M := 883 nb := 64 conf := gomas.NewConf() conf.LB = 0 zeromean := cmat.NewFloatNormSource() A := cmat.NewMatrix(M, N) A.SetFrom(zeromean) tauq := cmat.NewMatrix(N, 1) taup := cmat.NewMatrix(N, 1) A1 := cmat.NewCopy(A) tauq1 := cmat.NewMatrix(N, 1) taup1 := cmat.NewMatrix(N, 1) W := lapackd.Workspace(M + N) W1 := lapackd.Workspace(nb * (M + N + 1)) lapackd.BDReduce(A, tauq, taup, W, conf) conf.LB = nb lapackd.BDReduce(A1, tauq1, taup1, W1, conf) // unblk.BiRed(A) == blk.BiRed(A) blasd.Plus(A1, A, 1.0, -1.0, gomas.NONE) blasd.Axpy(tauq1, tauq, -1.0) blasd.Axpy(taup1, taup, -1.0) nrm := lapackd.NormP(A1, lapackd.NORM_ONE) t.Logf("M=%d, N=%d || unblk.BiRed(A) - blk.BiRed(A)||_1 : %e\n", M, N, nrm) nrm = lapackd.NormP(taup1, lapackd.NORM_ONE) t.Logf(" || unblk.BiRed(A).tauq - blk.BiRed(A).taup||_1 : %e\n", nrm) nrm = lapackd.NormP(tauq1, lapackd.NORM_ONE) t.Logf(" || unblk.BiRed(A).taup - blk.BiRed(A).tauq||_1 : %e\n", nrm) }
func test_trdevd(N, flags, kind int, verbose bool, t *testing.T) { var At, sD, sE, tmp cmat.FloatMatrix A0 := cmat.NewMatrix(N, N) desc := setTrdDiagonals(A0, kind) At.SubMatrix(A0, 0, 0, N, N) sD.Diag(A0, 0) sE.Diag(A0, 1) D := cmat.NewCopy(&sD) E := cmat.NewCopy(&sE) V := cmat.NewMatrix(N, N) sD.Diag(V, 0) sD.Add(1.0) W := cmat.NewMatrix(4*N, 1) C := cmat.NewMatrix(N, N) if verbose && N < 10 { t.Logf("A0:\n%v\n", A0.ToString("%6.3f")) t.Logf("V.pre:\n%v\n", V.ToString("%6.3f")) } lapackd.TRDEigen(D, E, V, W, flags|gomas.WANTV) for k := 0; k < N-1; k++ { if E.GetAt(k) != 0.0 { t.Logf("E[%d] != 0.0 (%e)\n", k, E.GetAt(k)) } } blasd.Mult(C, V, V, 1.0, 0.0, gomas.TRANSB) sD.Diag(C) sD.Add(-1.0) nrmv := lapackd.NormP(C, lapackd.NORM_ONE) blasd.Mult(C, V, A0, 1.0, 0.0, gomas.TRANSA) blasd.Mult(&At, C, V, 1.0, 0.0, gomas.NONE) if verbose && N < 10 { t.Logf("D:\n%v\n", asRow(&tmp, D).ToString("%6.3f")) t.Logf("V:\n%v\n", V.ToString("%6.3f")) t.Logf("V.T*A*V\n%v\n", At.ToString("%6.3f")) } sD.Diag(&At) blasd.Axpy(&sD, D, -1.0) nrma := lapackd.NormP(&At, lapackd.NORM_ONE) t.Logf("N=%d [%s] ||V.T*A*V - eigen(A)||_1: %e\n", N, desc, nrma) t.Logf(" ||I - V.T*V||_1: %e\n", nrmv) }
// test: unblk.QLFactor == blk.QLFactor func TestQLFactor(t *testing.T) { var t0 cmat.FloatMatrix M := 911 N := 835 nb := 32 conf := gomas.NewConf() A := cmat.NewMatrix(M, N) src := cmat.NewFloatNormSource() A.SetFrom(src) tau := cmat.NewMatrix(M, 1) W := cmat.NewMatrix(M+N, 1) A1 := cmat.NewCopy(A) tau1 := cmat.NewCopy(tau) conf.LB = 0 lapackd.QLFactor(A, tau, W, conf) conf.LB = nb W1 := lapackd.Workspace(lapackd.QLFactorWork(A1, conf)) lapackd.QLFactor(A1, tau1, W1, conf) if N < 10 { t.Logf("unblkQL(A):\n%v\n", A) t0.SetBuf(1, tau.Len(), 1, tau.Data()) t.Logf("unblkQL.tau:\n%v\n", &t0) t.Logf("blkQL(A):\n%v\n", A1) t0.SetBuf(1, tau1.Len(), 1, tau1.Data()) t.Logf("blkQL.tau:\n%v\n", &t0) } blasd.Plus(A1, A, 1.0, -1.0, gomas.NONE) nrm := lapackd.NormP(A1, lapackd.NORM_ONE) t.Logf("M=%d, N=%d ||blkQL(A) - unblkQL(A)||_1: %e\n", M, N, nrm) blasd.Axpy(tau1, tau, -1.0) nrm = blasd.Nrm2(tau1) t.Logf(" ||blkQL.tau - unblkQL.tau||_1: %e\n", nrm) }
/* * Update vector with compact WY Householder block * (I - Y*T*Y.T)*v = v - Y*T*Y.T*v * * LEFT: * 1 | 0 * v0 = v0 = v0 * 0 | Q v1 Q*v1 = v1 - Y*T*Y.T*v1 * * 1 | 0 * v0 = v0 = v0 * 0 | Q.T v1 Q.T*v1 = v1 - Y*T.T*Y.T*v1 * * RIGHT: * v0 | v1 * 1 | 0 = v0 | v1*Q = v0 | v1 - v1*Y*T*Y.T * 0 | Q * * v0 | v1 * 1 | 0 = v0 | v1*Q.T = v0 | v1 - v1*Y*T.T*Y.T * 0 | Q.T */ func updateVecLeftWY2(v, Y1, Y2, T, w *cmat.FloatMatrix, bits int) { var v1, v2 cmat.FloatMatrix var w0 cmat.FloatMatrix v1.SubMatrix(v, 1, 0, n(Y1), 1) v2.SubMatrix(v, n(Y1)+1, 0, m(Y2), 1) w0.SubMatrix(w, 0, 0, m(Y1), 1) // w0 := Y1.T*v1 + Y2.T*v2 blasd.Copy(&w0, &v1) blasd.MVMultTrm(&w0, Y1, 1.0, gomas.LOWER|gomas.UNIT|gomas.TRANS) blasd.MVMult(&w0, Y2, &v2, 1.0, 1.0, gomas.TRANS) // w0 := op(T)*w0 blasd.MVMultTrm(&w0, T, 1.0, bits|gomas.UPPER) // v2 := v2 - Y2*w0 blasd.MVMult(&v2, Y2, &w0, -1.0, 1.0, gomas.NONE) // v1 := v1 - Y1*w0 blasd.MVMultTrm(&w0, Y1, 1.0, gomas.LOWER|gomas.UNIT) blasd.Axpy(&v1, &w0, -1.0) }
/* From LAPACK/dlarf.f * * Applies a real elementary reflector H to a real m by n matrix A, * from either the left or the right. H is represented in the form * * H = I - tau * ( 1 ) * ( 1 v.T ) * ( v ) * * where tau is a real scalar and v is a real vector. * * If tau = 0, then H is taken to be the unit cmat. * * A is /a1\ a1 := a1 - w1 * \A2/ A2 := A2 - v*w1 * w1 := tau*(a1 + A2.T*v) if side == LEFT * := tau*(a1 + A2*v) if side == RIGHT * * Allocates/frees intermediate work space matrix w1. */ func applyHouseholder(tau, v, a1, A2 *cmat.FloatMatrix, flags int) { tval := tau.Get(0, 0) if tval == 0.0 { return } w1 := cmat.NewCopy(a1) if flags&gomas.LEFT != 0 { // w1 = a1 + A2.T*v blasd.MVMult(w1, A2, v, 1.0, 1.0, gomas.TRANSA) } else { // w1 = a1 + A2*v blasd.MVMult(w1, A2, v, 1.0, 1.0, gomas.NONE) } // w1 = tau*w1 blasd.Scale(w1, tval) // a1 = a1 - w1 blasd.Axpy(a1, w1, -1.0) // A2 = A2 - v*w1 blasd.MVUpdate(A2, v, w1, -1.0) }
/* * This is adaptation of TRIRED_LAZY_UNB algorithm from (1). */ func unblkBuildTridiagUpper(A, tauq, Y, W *cmat.FloatMatrix) { var ATL, ABR cmat.FloatMatrix var A00, a01, A02, a11, a12, A22 cmat.FloatMatrix var YTL, YBR cmat.FloatMatrix var Y00, y01, Y02, y11, y12, Y22 cmat.FloatMatrix var tqT, tqB, tq0, tauq1, tq2 cmat.FloatMatrix var w12 cmat.FloatMatrix var v0 float64 util.Partition2x2( &ATL, nil, nil, &ABR, A, 0, 0, util.PBOTTOMRIGHT) util.Partition2x2( &YTL, nil, nil, &YBR, Y, 0, 0, util.PBOTTOMRIGHT) util.Partition2x1( &tqT, &tqB, tauq, 0, util.PBOTTOM) k := 0 for k < n(Y) { util.Repartition2x2to3x3(&ATL, &A00, &a01, &A02, nil, &a11, &a12, nil, nil, &A22, A, 1, util.PTOPLEFT) util.Repartition2x2to3x3(&YTL, &Y00, &y01, &Y02, nil, &y11, &y12, nil, nil, &Y22, Y, 1, util.PTOPLEFT) util.Repartition2x1to3x1(&tqT, &tq0, &tauq1, &tq2, tauq, 1, util.PTOP) // set temp vectors for this round w12.SubMatrix(Y, -1, 0, 1, n(&Y02)) // ------------------------------------------------------ if n(&Y02) > 0 { aa := blasd.Dot(&a12, &y12) aa += blasd.Dot(&y12, &a12) a11.Set(0, 0, a11.Get(0, 0)-aa) // a01 := a01 - A02*y12 blasd.MVMult(&a01, &A02, &y12, -1.0, 1.0, gomas.NONE) // a01 := a01 - Y02*a12 blasd.MVMult(&a01, &Y02, &a12, -1.0, 1.0, gomas.NONE) // restore superdiagonal value a12.Set(0, 0, v0) } // Compute householder to zero subdiagonal entries computeHouseholderRev(&a01, &tauq1) tauqv := tauq1.Get(0, 0) // set sub&iagonal to unit v0 = a01.Get(-1, 0) a01.Set(-1, 0, 1.0) // y01 := tauq*A00*a01 blasd.MVMultSym(&y01, &A00, &a01, tauqv, 0.0, gomas.UPPER) // w12 := A02.T*a01 blasd.MVMult(&w12, &A02, &a01, 1.0, 0.0, gomas.TRANS) // y01 := y01 - Y02*(A02.T*a01) blasd.MVMult(&y01, &Y02, &w12, -tauqv, 1.0, gomas.NONE) // w12 := Y02.T*a01 blasd.MVMult(&w12, &Y02, &a01, 1.0, 0.0, gomas.TRANS) // y01 := y01 - A02*(Y02.T*a01) blasd.MVMult(&y01, &A02, &w12, -tauqv, 1.0, gomas.NONE) // beta := tauq*a01.T*y01 beta := tauqv * blasd.Dot(&a01, &y01) // y01 := y01 - 0.5*beta*a01 blasd.Axpy(&y01, &a01, -0.5*beta) // ------------------------------------------------------ k += 1 util.Continue3x3to2x2( &ATL, nil, nil, &ABR, &A00, &a11, &A22, A, util.PTOPLEFT) util.Continue3x3to2x2( &YTL, nil, nil, &YBR, &Y00, &y11, &Y22, A, util.PTOPLEFT) util.Continue3x1to2x1( &tqT, &tqB, &tq0, &tauq1, tauq, util.PTOP) } // restore superdiagonal value A.Set(m(&ATL)-1, n(&ATL), v0) }
/* * This is adaptation of TRIRED_LAZY_UNB algorithm from (1). */ func unblkBuildTridiagLower(A, tauq, Y, W *cmat.FloatMatrix) { var ATL, ABR cmat.FloatMatrix var A00, a10, a11, A20, a21, A22 cmat.FloatMatrix var YTL, YBR cmat.FloatMatrix var Y00, y10, y11, Y20, y21, Y22 cmat.FloatMatrix var tqT, tqB, tq0, tauq1, tq2 cmat.FloatMatrix var w12 cmat.FloatMatrix var v0 float64 util.Partition2x2( &ATL, nil, nil, &ABR, A, 0, 0, util.PTOPLEFT) util.Partition2x2( &YTL, nil, nil, &YBR, Y, 0, 0, util.PTOPLEFT) util.Partition2x1( &tqT, &tqB, tauq, 0, util.PTOP) k := 0 for k < n(Y) { util.Repartition2x2to3x3(&ATL, &A00, nil, nil, &a10, &a11, nil, &A20, &a21, &A22, A, 1, util.PBOTTOMRIGHT) util.Repartition2x2to3x3(&YTL, &Y00, nil, nil, &y10, &y11, nil, &Y20, &y21, &Y22, Y, 1, util.PBOTTOMRIGHT) util.Repartition2x1to3x1(&tqT, &tq0, &tauq1, &tq2, tauq, 1, util.PBOTTOM) // set temp vectors for this round //w12.SetBuf(y10.Len(), 1, y10.Len(), W.Data()) w12.SubMatrix(Y, 0, 0, 1, n(&Y00)) // ------------------------------------------------------ if n(&Y00) > 0 { aa := blasd.Dot(&a10, &y10) aa += blasd.Dot(&y10, &a10) a11.Set(0, 0, a11.Get(0, 0)-aa) // a21 := a21 - A20*y10 blasd.MVMult(&a21, &A20, &y10, -1.0, 1.0, gomas.NONE) // a21 := a21 - Y20*a10 blasd.MVMult(&a21, &Y20, &a10, -1.0, 1.0, gomas.NONE) // restore subdiagonal value a10.Set(0, -1, v0) } // Compute householder to zero subdiagonal entries computeHouseholderVec(&a21, &tauq1) tauqv := tauq1.Get(0, 0) // set subdiagonal to unit v0 = a21.Get(0, 0) a21.Set(0, 0, 1.0) // y21 := tauq*A22*a21 blasd.MVMultSym(&y21, &A22, &a21, tauqv, 0.0, gomas.LOWER) // w12 := A20.T*a21 blasd.MVMult(&w12, &A20, &a21, 1.0, 0.0, gomas.TRANS) // y21 := y21 - Y20*(A20.T*a21) blasd.MVMult(&y21, &Y20, &w12, -tauqv, 1.0, gomas.NONE) // w12 := Y20.T*a21 blasd.MVMult(&w12, &Y20, &a21, 1.0, 0.0, gomas.TRANS) // y21 := y21 - A20*(Y20.T*a21) blasd.MVMult(&y21, &A20, &w12, -tauqv, 1.0, gomas.NONE) // beta := tauq*a21.T*y21 beta := tauqv * blasd.Dot(&a21, &y21) // y21 := y21 - 0.5*beta*a21 blasd.Axpy(&y21, &a21, -0.5*beta) // ------------------------------------------------------ k += 1 util.Continue3x3to2x2( &ATL, nil, nil, &ABR, &A00, &a11, &A22, A, util.PBOTTOMRIGHT) util.Continue3x3to2x2( &YTL, nil, nil, &YBR, &Y00, &y11, &Y22, A, util.PBOTTOMRIGHT) util.Continue3x1to2x1( &tqT, &tqB, &tq0, &tauq1, tauq, util.PBOTTOM) } // restore subdiagonal value A.Set(m(&ATL), n(&ATL)-1, v0) }
/* * Compute unblocked bidiagonal reduction for A when M >= N * * Diagonal and first super/sub diagonal are overwritten with the * upper/lower bidiagonal matrix B. * * This computing (1-tauq*v*v.T)*A*(1-taup*u.u.T) from left to right. */ func unblkReduceBidiagLeft(A, tauq, taup, W *cmat.FloatMatrix) { var ATL, ABR cmat.FloatMatrix var A00, a11, a12t, a21, A22 cmat.FloatMatrix var tqT, tqB, tq0, tauq1, tq2 cmat.FloatMatrix var tpT, tpB, tp0, taup1, tp2 cmat.FloatMatrix var y21, z21 cmat.FloatMatrix var v0 float64 util.Partition2x2( &ATL, nil, nil, &ABR, A, 0, 0, util.PTOPLEFT) util.Partition2x1( &tqT, &tqB, tauq, 0, util.PTOP) util.Partition2x1( &tpT, &tpB, taup, 0, util.PTOP) for m(&ABR) > 0 && n(&ABR) > 0 { util.Repartition2x2to3x3(&ATL, &A00, nil, nil, nil, &a11, &a12t, nil, &a21, &A22, A, 1, util.PBOTTOMRIGHT) util.Repartition2x1to3x1(&tqT, &tq0, &tauq1, &tq2, tauq, 1, util.PBOTTOM) util.Repartition2x1to3x1(&tpT, &tp0, &taup1, &tp2, taup, 1, util.PBOTTOM) // set temp vectors for this round y21.SetBuf(n(&a12t), 1, n(&a12t), W.Data()) z21.SetBuf(m(&a21), 1, m(&a21), W.Data()[y21.Len():]) // ------------------------------------------------------ // Compute householder to zero subdiagonal entries computeHouseholder(&a11, &a21, &tauq1) // y21 := a12 + A22.T*a21 blasd.Axpby(&y21, &a12t, 1.0, 0.0) blasd.MVMult(&y21, &A22, &a21, 1.0, 1.0, gomas.TRANSA) // a12t := a12t - tauq*y21 tauqv := tauq1.Get(0, 0) blasd.Axpy(&a12t, &y21, -tauqv) // Compute householder to zero elements above 1st superdiagonal computeHouseholderVec(&a12t, &taup1) v0 = a12t.Get(0, 0) a12t.Set(0, 0, 1.0) taupv := taup1.Get(0, 0) // [v == a12t, u == a21] beta := blasd.Dot(&y21, &a12t) // z21 := tauq*beta*u blasd.Axpby(&z21, &a21, tauqv*beta, 0.0) // z21 := A22*v - z21 blasd.MVMult(&z21, &A22, &a12t, 1.0, -1.0, gomas.NONE) // A22 := A22 - tauq*u*y21 blasd.MVUpdate(&A22, &a21, &y21, -tauqv) // A22 := A22 - taup*z21*v blasd.MVUpdate(&A22, &z21, &a12t, -taupv) a12t.Set(0, 0, v0) // ------------------------------------------------------ util.Continue3x3to2x2( &ATL, nil, nil, &ABR, &A00, &a11, &A22, A, util.PBOTTOMRIGHT) util.Continue3x1to2x1( &tqT, &tqB, &tq0, &tauq1, tauq, util.PBOTTOM) util.Continue3x1to2x1( &tpT, &tpB, &tp0, &taup1, taup, util.PBOTTOM) } }
/* * This is adaptation of BIRED_LAZY_UNB algorithm from (1). * * Z matrix accumulates updates of row transformations i.e. first * Householder that zeros off diagonal entries on row. Vector z21 * is updates for current round, Z20 are already accumulated updates. * Vector z21 updates a12 before next transformation. * * Y matrix accumulates updates on column tranformations ie Householder * that zeros elements below sub-diagonal. Vector y21 is updates for current * round, Y20 are already accumulated updates. Vector y21 updates * a21 befor next transformation. * * Z, Y matrices upper trigonal part is not needed, temporary vector * w00 that has maximum length of n(Y) is placed on the last column of * Z matrix on each iteration. */ func unblkBuildBidiagRight(A, tauq, taup, Y, Z *cmat.FloatMatrix) { var ATL, ABL, ABR cmat.FloatMatrix var A00, a01, A02, a10, a11, a12t, A20, a21, A22 cmat.FloatMatrix var YTL, YBR, ZTL, ZBR cmat.FloatMatrix var Y00, y10, Y20, y11, y21, Y22 cmat.FloatMatrix var Z00, z10, Z20, z11, z21, Z22 cmat.FloatMatrix var tqT, tqB, tq0, tauq1, tq2 cmat.FloatMatrix var tpT, tpB, tp0, taup1, tp2 cmat.FloatMatrix var w00 cmat.FloatMatrix var v0 float64 // Y is workspace for building updates for first Householder. // And Z is space for build updates for second Householder // Y is n(A)-2,nb and Z is m(A)-1,nb util.Partition2x2( &ATL, nil, &ABL, &ABR, A, 0, 0, util.PTOPLEFT) util.Partition2x2( &YTL, nil, nil, &YBR, Y, 0, 0, util.PTOPLEFT) util.Partition2x2( &ZTL, nil, nil, &ZBR, Z, 0, 0, util.PTOPLEFT) util.Partition2x1( &tqT, &tqB, tauq, 0, util.PTOP) util.Partition2x1( &tpT, &tpB, taup, 0, util.PTOP) k := 0 for k < n(Y) { util.Repartition2x2to3x3(&ATL, &A00, &a01, &A02, &a10, &a11, &a12t, &A20, &a21, &A22, A, 1, util.PBOTTOMRIGHT) util.Repartition2x2to3x3(&YTL, &Y00, nil, nil, &y10, &y11, nil, &Y20, &y21, &Y22, Y, 1, util.PBOTTOMRIGHT) util.Repartition2x2to3x3(&ZTL, &Z00, nil, nil, &z10, &z11, nil, &Z20, &z21, &Z22, Z, 1, util.PBOTTOMRIGHT) util.Repartition2x1to3x1(&tqT, &tq0, &tauq1, &tq2, tauq, 1, util.PBOTTOM) util.Repartition2x1to3x1(&tpT, &tp0, &taup1, &tp2, taup, 1, util.PBOTTOM) // set temp vectors for this round, w00.SubMatrix(Z, 0, n(Z)-1, m(&A02), 1) // ------------------------------------------------------ // u10 == a10, U20 == A20, u21 == a21, // v10 == a01, V20 == A02, v21 == a12t if n(&Y20) > 0 { // a11 := a11 - u10t*z10 - y10*v10 aa := blasd.Dot(&a10, &z10) aa += blasd.Dot(&y10, &a01) a11.Set(0, 0, a11.Get(0, 0)-aa) // a12t := a12t - V20*z10 - Z20*u10 blasd.MVMult(&a12t, &A02, &y10, -1.0, 1.0, gomas.TRANS) blasd.MVMult(&a12t, &Z20, &a10, -1.0, 1.0, gomas.NONE) // a21 := a21 - Y20*v10 - U20*z10 blasd.MVMult(&a21, &Y20, &a01, -1.0, 1.0, gomas.NONE) blasd.MVMult(&a21, &A20, &z10, -1.0, 1.0, gomas.NONE) // here restore bidiagonal entry a10.Set(0, -1, v0) } // Compute householder to zero superdiagonal entries computeHouseholder(&a11, &a12t, &taup1) taupv := taup1.Get(0, 0) // y21 := a21 + A22*v21 - Y20*U20.T*v21 - V20*Z20.T*v21 blasd.Axpby(&y21, &a21, 1.0, 0.0) blasd.MVMult(&y21, &A22, &a12t, 1.0, 1.0, gomas.NONE) // w00 := U20.T*v21 [= A02*a12t] blasd.MVMult(&w00, &A02, &a12t, 1.0, 0.0, gomas.NONE) // y21 := y21 - U20*w00 [U20 == A20] blasd.MVMult(&y21, &Y20, &w00, -1.0, 1.0, gomas.NONE) // w00 := Z20.T*v21 blasd.MVMult(&w00, &Z20, &a12t, 1.0, 0.0, gomas.TRANS) // y21 := y21 - V20*w00 [V20 == A02.T] blasd.MVMult(&y21, &A20, &w00, -1.0, 1.0, gomas.NONE) // a21 := a21 - taup*y21 blasd.Scale(&y21, taupv) blasd.Axpy(&a21, &y21, -1.0) // Compute householder to zero elements below 1st subdiagonal computeHouseholderVec(&a21, &tauq1) v0 = a21.Get(0, 0) a21.Set(0, 0, 1.0) tauqv := tauq1.Get(0, 0) // z21 := tauq*(A22*y - V20*Y20.T*u - Z20*U20.T*u - beta*v) // [v == a12t, u == a21] beta := blasd.Dot(&y21, &a21) // z21 := beta*v blasd.Axpby(&z21, &a12t, beta, 0.0) // w00 = Y20.T*u blasd.MVMult(&w00, &Y20, &a21, 1.0, 0.0, gomas.TRANS) // z21 = z21 + V20*w00 == A02.T*w00 blasd.MVMult(&z21, &A02, &w00, 1.0, 1.0, gomas.TRANS) // w00 := U20.T*u (U20.T == A20.T) blasd.MVMult(&w00, &A20, &a21, 1.0, 0.0, gomas.TRANS) // z21 := z21 + Z20*w00 blasd.MVMult(&z21, &Z20, &w00, 1.0, 1.0, gomas.NONE) // z21 := -tauq*z21 + tauq*A22*v blasd.MVMult(&z21, &A22, &a21, tauqv, -tauqv, gomas.TRANS) // ------------------------------------------------------ k += 1 util.Continue3x3to2x2( &ATL, nil, &ABL, &ABR, &A00, &a11, &A22, A, util.PBOTTOMRIGHT) util.Continue3x3to2x2( &YTL, nil, nil, &YBR, &Y00, &y11, &Y22, Y, util.PBOTTOMRIGHT) util.Continue3x3to2x2( &ZTL, nil, nil, &ZBR, &Z00, &z11, &Z22, Z, util.PBOTTOMRIGHT) util.Continue3x1to2x1( &tqT, &tqB, &tq0, &tauq1, tauq, util.PBOTTOM) util.Continue3x1to2x1( &tpT, &tpB, &tp0, &taup1, taup, util.PBOTTOM) } // restore ABL.Set(0, -1, v0) }