Пример #1
0
/*
 * Multiply and replace C with Q*C or Q.T*C where Q is a real orthogonal matrix
 * defined as the product of k elementary reflectors and block reflector T
 *
 *    Q = H(1) H(2) . . . H(k)
 *
 * as returned by DecomposeQRT().
 *
 * Arguments:
 *  C     On entry, the M-by-N matrix C. On exit C is overwritten by Q*C or Q.T*C.
 *
 *  A     QR factorization as returned by QRTFactor() where the lower trapezoidal
 *        part holds the elementary reflectors.
 *
 *  T     The block reflector computed from elementary reflectors as returned by
 *        DecomposeQRT() or computed from elementary reflectors and scalar coefficients
 *        by BuildT()
 *
 *  W     Workspace, size as returned by QRTMultWork()
 *
 *  conf  Blocking configuration
 *
 *  flags Indicators. Valid indicators LEFT, RIGHT, TRANS, NOTRANS
 *
 * Preconditions:
 *   a.   cols(A) == cols(T),
 *          columns A define number of elementary reflector, must match order of block reflector.
 *   b.   if conf.LB == 0, cols(T) == rows(T)
 *          unblocked invocation, block reflector T is upper triangular
 *   c.   if conf.LB != 0, rows(T) == conf.LB
 *          blocked invocation, T is sequence of triangular block reflectors of order LB
 *   d.   if LEFT, rows(C) >= cols(A) && cols(C) >= rows(A)
 *
 *   e.   if RIGHT, cols(C) >= cols(A) && rows(C) >= rows(A)
 *
 * Compatible with lapack.DGEMQRT
 */
func QRTMult(C, A, T, W *cmat.FloatMatrix, flags int, confs ...*gomas.Config) *gomas.Error {
	var err *gomas.Error = nil
	conf := gomas.CurrentConf(confs...)

	wsz := QRTMultWork(C, T, flags, conf)
	if W == nil || W.Len() < wsz {
		return gomas.NewError(gomas.EWORK, "QRTMult", wsz)
	}
	ok := false
	switch flags & gomas.RIGHT {
	case gomas.RIGHT:
		ok = n(C) >= m(A)
	default:
		ok = m(C) >= n(A)
	}
	if !ok {
		return gomas.NewError(gomas.ESIZE, "QRTMult")
	}

	var Wrk cmat.FloatMatrix
	if flags&gomas.RIGHT != 0 {
		Wrk.SetBuf(m(C), conf.LB, m(C), W.Data())
		blockedMultQTRight(C, A, T, &Wrk, flags, conf)

	} else {
		Wrk.SetBuf(n(C), conf.LB, n(C), W.Data())
		blockedMultQTLeft(C, A, T, &Wrk, flags, conf)
	}
	return err
}
Пример #2
0
/*
 * Compute
 *   B = B*diag(D).-1      flags & RIGHT == true
 *   B = diag(D).-1*B      flags & LEFT  == true
 *
 * If flags is LEFT (RIGHT) then element-wise divides columns (rows) of B with vector D.
 *
 * Arguments:
 *   B     M-by-N matrix if flags&RIGHT == true or N-by-M matrix if flags&LEFT == true
 *
 *   D     N element column or row vector or N-by-N matrix
 *
 *   flags Indicator bits, LEFT or RIGHT
 */
func SolveDiag(B, D *cmat.FloatMatrix, flags int, confs ...*gomas.Config) *gomas.Error {
	var c, d0 cmat.FloatMatrix
	var d *cmat.FloatMatrix

	conf := gomas.CurrentConf(confs...)
	d = D
	if !D.IsVector() {
		d0.Diag(D)
		d = &d0
	}
	dn := d0.Len()
	br, bc := B.Size()
	switch flags & (gomas.LEFT | gomas.RIGHT) {
	case gomas.LEFT:
		if br != dn {
			return gomas.NewError(gomas.ESIZE, "SolveDiag")
		}
		// scale rows;
		for k := 0; k < dn; k++ {
			c.Row(B, k)
			blasd.InvScale(&c, d.GetAt(k), conf)
		}
	case gomas.RIGHT:
		if bc != dn {
			return gomas.NewError(gomas.ESIZE, "SolveDiag")
		}
		// scale columns
		for k := 0; k < dn; k++ {
			c.Column(B, k)
			blasd.InvScale(&c, d.GetAt(k), conf)
		}
	}
	return nil
}
Пример #3
0
// Compute i'th updated element of rank-1 update vector
func trdsecUpdateElemDelta(d, delta *cmat.FloatMatrix, index int, rho float64) float64 {
	var n0, n1, dn, dk, val, p0, p1 float64
	var k, N int

	N = d.Len()
	dk = d.GetAt(index)
	dn = delta.GetAt(N - 1)

	// compute; prod j; (lambda_j - d_k)/(d_j - d_k), j = 0 .. index-1
	p0 = 1.0
	for k = 0; k < index; k++ {
		n0 = delta.GetAt(k)
		n1 = d.GetAt(k) - dk
		p0 = p0 * (n0 / n1)
	}

	p1 = 1.0
	for k = index; k < N-1; k++ {
		n0 = delta.GetAt(k)
		n1 = d.GetAt(k+1) - dk
		p1 = p1 * (n0 / n1)
	}
	val = p0 * p1 * (dn / rho)
	return math.Sqrt(math.Abs(val))
}
Пример #4
0
func sortEigenVec(D, U, V, C *cmat.FloatMatrix, updown int) {
	var sD, m0, m1 cmat.FloatMatrix

	N := D.Len()
	for k := 0; k < N-1; k++ {
		sD.SubVector(D, k, N-k)
		pk := vecMinMax(&sD, -updown)
		if pk != 0 {
			t0 := D.GetAt(k)
			D.SetAt(k, D.GetAt(pk+k))
			D.SetAt(k+pk, t0)
			if U != nil {
				m0.Column(U, k)
				m1.Column(U, k+pk)
				blasd.Swap(&m1, &m0)
			}
			if V != nil {
				m0.Row(V, k)
				m1.Row(V, k+pk)
				blasd.Swap(&m1, &m0)
			}
			if C != nil {
				m0.Column(C, k)
				m1.Column(C, k+pk)
				blasd.Swap(&m1, &m0)
			}
		}
	}
}
Пример #5
0
/*
 * Compute QR factorization of a M-by-N matrix A using compact WY transformation: A = Q * R,
 * where Q = I - Y*T*Y.T, T is block reflector and Y holds elementary reflectors as lower
 * trapezoidal matrix saved below diagonal elements of the matrix A.
 *
 * Arguments:
 *  A    On entry, the M-by-N matrix A. On exit, the elements on and above
 *       the diagonal contain the min(M,N)-by-N upper trapezoidal matrix R.
 *       The elements below the diagonal with the matrix 'T', represent
 *       the ortogonal matrix Q as product of elementary reflectors.
 *
 * T     On exit, the K block reflectors which, together with trilu(A) represent
 *       the ortogonal matrix Q as Q = I - Y*T*Y.T where Y = trilu(A).
 *       K is ceiling(N/LB) where LB is blocking size from used blocking configuration.
 *       The matrix T is LB*N augmented matrix of K block reflectors,
 *       T = [T(0) T(1) .. T(K-1)].  Block reflector T(n) is LB*LB matrix, expect
 *       reflector T(K-1) that is IB*IB matrix  where IB = min(LB, K % LB)
 *
 * W     Workspace, required size returned by QRTFactorWork().
 *
 * conf  Optional blocking configuration. If not provided then default configuration
 *       is used.
 *
 * Returns:
 *      Error indicator.
 *
 * QRTFactor is compatible with lapack.DGEQRT
 */
func QRTFactor(A, T, W *cmat.FloatMatrix, confs ...*gomas.Config) *gomas.Error {
	var err *gomas.Error = nil
	conf := gomas.CurrentConf(confs...)
	ok := false
	rsize := 0

	if m(A) < n(A) {
		return gomas.NewError(gomas.ESIZE, "QRTFactor")
	}
	wsz := QRTFactorWork(A, conf)
	if W == nil || W.Len() < wsz {
		return gomas.NewError(gomas.EWORK, "QRTFactor", wsz)
	}

	tr, tc := T.Size()
	if conf.LB == 0 || conf.LB > n(A) {
		ok = tr == tc && tr == n(A)
		rsize = n(A) * n(A)
	} else {
		ok = tr == conf.LB && tc == n(A)
		rsize = conf.LB * n(A)
	}
	if !ok {
		return gomas.NewError(gomas.ESMALL, "QRTFactor", rsize)
	}

	if conf.LB == 0 || n(A) <= conf.LB {
		err = unblockedQRT(A, T, W)
	} else {
		Wrk := cmat.MakeMatrix(n(A), conf.LB, W.Data())
		err = blockedQRT(A, T, Wrk, conf)
	}
	return err
}
Пример #6
0
/*
 * Compute LDL^T factorization of real symmetric matrix.
 *
 * Computes of a real symmetric matrix A using Bunch-Kauffman pivoting method.
 * The form of factorization is
 *
 *    A = L*D*L.T  or A = U*D*U.T
 *
 * where L (or U) is product of permutation and unit lower (or upper) triangular matrix
 * and D is block diagonal symmetric matrix with 1x1 and 2x2 blocks.
 *
 * Arguments
 *  A     On entry, the N-by-N symmetric matrix A. If flags bit LOWER (or UPPER) is set then
 *        lower (or upper) triangular matrix and strictly upper (or lower) part is not
 *        accessed. On exit, the block diagonal matrix D and lower (or upper) triangular
 *        product matrix L (or U).
 *
 *  W     Workspace, size as returned by WorksizeBK().
 *
 *  ipiv  Pivot vector. On exit details of interchanges and the block structure of D. If
 *        ipiv[k] > 0 then D[k,k] is 1x1 and rows and columns k and ipiv[k]-1 were changed.
 *        If ipiv[k] == ipiv[k+1] < 0 then D[k,k] is 2x2. If A is lower then rows and
 *        columns k+1 and ipiv[k]-1  were changed. And if A is upper then rows and columns
 *        k and ipvk[k]-1 were changed.
 *
 *  flags Indicator bits, LOWER or UPPER.
 *
 *  confs Optional blocking configuration. If not provided then default blocking
 *        as returned by DefaultConf() is used.
 *
 *  Unblocked algorithm is used if blocking configuration LB is zero or if N < LB.
 *
 *  Compatible with lapack.SYTRF.
 */
func BKFactor(A, W *cmat.FloatMatrix, ipiv Pivots, flags int, confs ...*gomas.Config) *gomas.Error {
	var err *gomas.Error = nil
	conf := gomas.CurrentConf(confs...)

	for k, _ := range ipiv {
		ipiv[k] = 0
	}
	wsz := BKFactorWork(A, conf)
	if W.Len() < wsz {
		return gomas.NewError(gomas.EWORK, "DecomposeBK", wsz)
	}

	var Wrk cmat.FloatMatrix
	if n(A) < conf.LB || conf.LB == 0 {
		// make workspace rows(A)*2 matrix
		Wrk.SetBuf(m(A), 2, m(A), W.Data())
		if flags&gomas.LOWER != 0 {
			err, _ = unblkDecompBKLower(A, &Wrk, ipiv, conf)
		} else if flags&gomas.UPPER != 0 {
			err, _ = unblkDecompBKUpper(A, &Wrk, ipiv, conf)
		}
	} else {
		// make workspace rows(A)*(LB+1) matrix
		Wrk.SetBuf(m(A), conf.LB+1, m(A), W.Data())
		if flags&gomas.LOWER != 0 {
			err = blkDecompBKLower(A, &Wrk, &ipiv, conf)
		} else if flags&gomas.UPPER != 0 {
			err = blkDecompBKUpper(A, &Wrk, &ipiv, conf)
		}
	}
	return err
}
Пример #7
0
/*
 * Reduce a general M-by-N matrix A to upper or lower bidiagonal form B
 * by an ortogonal transformation A = Q*B*P.T,  B = Q.T*A*P
 *
 *
 * Arguments
 *   A     On entry, the real M-by-N matrix. On exit the upper/lower
 *         bidiagonal matrix and ortogonal matrices Q and P.
 *
 *   tauq  Scalar factors for elementary reflector forming the
 *         ortogonal matrix Q.
 *
 *   taup  Scalar factors for elementary reflector forming the
 *         ortogonal matrix P.
 *
 *   W     Workspace needed for reduction.
 *
 *   conf  Current blocking configuration. Optional.
 *
 *
 * Details
 *
 * Matrices Q and P are products of elementary reflectors H(k) and G(k)
 *
 * If M > N:
 *     Q = H(1)*H(2)*...*H(N)   and P = G(1)*G(2)*...*G(N-1)
 *
 * where H(k) = 1 - tauq*u*u.T and G(k) = 1 - taup*v*v.T
 *
 * Elementary reflector H(k) are stored on columns of A below the diagonal with
 * implicit unit value on diagonal entry. Vector TAUQ holds corresponding scalar
 * factors. Reflector G(k) are stored on rows of A right of first superdiagonal
 * with implicit unit value on superdiagonal. Corresponding scalar factors are
 * stored on vector TAUP.
 *
 * If M < N:
 *   Q = H(1)*H(2)*...*H(N-1)   and P = G(1)*G(2)*...*G(N)
 *
 * where H(k) = 1 - tauq*u*u.T and G(k) = 1 - taup*v*v.T
 *
 * Elementary reflector H(k) are stored on columns of A below the first sub diagonal
 * with implicit unit value on sub diagonal entry. Vector TAUQ holds corresponding
 * scalar factors. Reflector G(k) are sotre on rows of A right of diagonal with
 * implicit unit value on superdiagonal. Corresponding scalar factors are stored
 * on vector TAUP.
 *
 * Contents of matrix A after reductions are as follows.
 *
 *    M = 6 and N = 5:                  M = 5 and N = 6:
 *
 *    (  d   e   v1  v1  v1 )           (  d   v1  v1  v1  v1  v1 )
 *    (  u1  d   e   v2  v2 )           (  e   d   v2  v2  v2  v2 )
 *    (  u1  u2  d   e   v3 )           (  u1  e   d   v3  v3  v3 )
 *    (  u1  u2  u3  d   e  )           (  u1  u2  e   d   v4  v4 )
 *    (  u1  u2  u3  u4  d  )           (  u1  u2  u3  e   d   v5 )
 *    (  u1  u2  u3  u4  u5 )
 */
func BDReduce(A, tauq, taup, W *cmat.FloatMatrix, confs ...*gomas.Config) *gomas.Error {
	var err *gomas.Error = nil
	conf := gomas.CurrentConf(confs...)
	_ = conf

	wmin := wsBired(A, 0)
	wsz := W.Len()
	if wsz < wmin {
		return gomas.NewError(gomas.EWORK, "ReduceBidiag", wmin)
	}
	lb := conf.LB
	wneed := wsBired(A, lb)
	if wneed > wsz {
		lb = estimateLB(A, wsz, wsBired)
	}
	if m(A) >= n(A) {
		if lb > 0 && n(A) > lb {
			blkBidiagLeft(A, tauq, taup, W, lb, conf)
		} else {
			unblkReduceBidiagLeft(A, tauq, taup, W)
		}
	} else {
		if lb > 0 && m(A) > lb {
			blkBidiagRight(A, tauq, taup, W, lb, conf)
		} else {
			unblkReduceBidiagRight(A, tauq, taup, W)
		}
	}
	return err
}
Пример #8
0
func MVMult(Y, A, X *cmat.FloatMatrix, alpha, beta float64, bits int, confs ...*gomas.Config) *gomas.Error {
	ok := true
	yr, yc := Y.Size()
	ar, ac := A.Size()
	xr, xc := X.Size()

	if ar*ac == 0 {
		return nil
	}
	if yr != 1 && yc != 1 {
		return gomas.NewError(gomas.ENEED_VECTOR, "MVMult")
	}
	if xr != 1 && xc != 1 {
		return gomas.NewError(gomas.ENEED_VECTOR, "MVMult")
	}
	nx := X.Len()
	ny := Y.Len()

	if bits&gomas.TRANSA != 0 {
		bits |= gomas.TRANS
	}
	if bits&gomas.TRANS != 0 {
		ok = ny == ac && nx == ar
	} else {
		ok = ny == ar && nx == ac
	}
	if !ok {
		return gomas.NewError(gomas.ESIZE, "MVMult")
	}
	if beta != 1.0 {
		vscal(Y, beta, ny)
	}
	gemv(Y, A, X, alpha, beta, bits, 0, nx, 0, ny)
	return nil
}
Пример #9
0
/*
 * Reduce general matrix A to upper Hessenberg form H by similiarity
 * transformation H = Q.T*A*Q.
 *
 * Arguments:
 *  A    On entry, the general matrix A. On exit, the elements on and
 *       above the first subdiagonal contain the reduced matrix H.
 *       The elements below the first subdiagonal with the vector tau
 *       represent the ortogonal matrix A as product of elementary reflectors.
 *
 *  tau  On exit, the scalar factors of the elementary reflectors.
 *
 *  W    Workspace, as defined by HessReduceWork()
 *
 *  conf The blocking configration.
 *
 * HessReduce is compatible with lapack.DGEHRD.
 */
func HessReduce(A, tau, W *cmat.FloatMatrix, confs ...*gomas.Config) *gomas.Error {
	var err *gomas.Error = nil
	conf := gomas.CurrentConf(confs...)

	wmin := m(A)
	wopt := HessReduceWork(A, conf)
	wsz := W.Len()
	if wsz < wmin {
		return gomas.NewError(gomas.EWORK, "ReduceHess", wmin)
	}
	// use blocked version if workspace big enough for blocksize 4
	lb := conf.LB
	if wsz < wopt {
		lb = estimateLB(A, wsz, wsHess)
	}
	if lb == 0 || n(A) <= lb {
		unblkHessGQvdG(A, tau, W, 0)
	} else {
		// blocked version
		var W0 cmat.FloatMatrix
		// shape workspace for blocked algorithm
		W0.SetBuf(m(A)+lb, lb, m(A)+lb, W.Data())
		blkHessGQvdG(A, tau, &W0, lb, conf)
	}
	return err
}
Пример #10
0
// d = |d| - |s|
func absMinus(d, s *cmat.FloatMatrix) *cmat.FloatMatrix {
	for k := 0; k < d.Len(); k++ {
		tmp := math.Abs(d.GetAt(k))
		d.SetAt(k, math.Abs(s.GetAt(k))-tmp)
	}
	return d
}
Пример #11
0
/*
 * Symmetric matrix-vector multiplication. Y = beta*Y + alpha*A*X
 */
func MVMultSym(Y, A, X *cmat.FloatMatrix, alpha, beta float64, bits int, confs ...*gomas.Config) *gomas.Error {
	ok := true
	yr, yc := Y.Size()
	ar, ac := A.Size()
	xr, xc := X.Size()

	if ar*ac == 0 {
		return nil
	}
	if yr != 1 && yc != 1 {
		return gomas.NewError(gomas.ENEED_VECTOR, "MVMultSym")
	}
	if xr != 1 && xc != 1 {
		return gomas.NewError(gomas.ENEED_VECTOR, "MVMultSym")
	}
	nx := X.Len()
	ny := Y.Len()

	ok = ny == ar && nx == ac && ac == ar
	if !ok {
		return gomas.NewError(gomas.ESIZE, "MVMultSym")
	}
	if beta != 1.0 {
		vscal(Y, beta, ny)
	}
	symv(Y, A, X, alpha, bits, nx)
	return nil
}
Пример #12
0
/*
 * Bidiagonal top to bottom implicit zero shift QR sweep
 */
func bdQRzero(D, E, Cr, Sr, Cl, Sl *cmat.FloatMatrix, saves bool) int {
	var d1, e1, d2, cosr, sinr, cosl, sinl, r float64
	N := D.Len()

	d1 = D.GetAtUnsafe(0)
	cosr = 1.0
	cosl = 1.0
	for k := 0; k < N-1; k++ {
		e1 = E.GetAtUnsafe(k)
		d2 = D.GetAtUnsafe(k + 1)
		cosr, sinr, r = ComputeGivens(d1*cosr, e1)
		if k > 0 {
			E.SetAtUnsafe(k-1, sinl*r)
		}
		cosl, sinl, r = ComputeGivens(cosl*r, sinr*d2)
		D.SetAtUnsafe(k, r)
		d1 = d2
		if saves {
			Cr.SetAtUnsafe(k, cosr)
			Sr.SetAtUnsafe(k, sinr)
			Cl.SetAtUnsafe(k, cosl)
			Sl.SetAtUnsafe(k, sinl)
		}
	}
	d2 = cosr * d2
	D.SetAtUnsafe(N-1, d2*cosl)
	E.SetAtUnsafe(N-2, d2*sinl)
	return N - 1
}
Пример #13
0
// Compute the updated rank-1 update vector with precomputed deltas
func trdsecUpdateVecDelta(z, Q, d *cmat.FloatMatrix, rho float64) {
	var delta cmat.FloatMatrix
	for i := 0; i < d.Len(); i++ {
		delta.Column(Q, i)
		zk := trdsecUpdateElemDelta(d, &delta, i, rho)
		z.SetAt(i, zk)
	}
}
Пример #14
0
// Compute eigenmatrix Q for updated eigenvalues in 'dl'.
func trdsecEigenBuild(Q, z, Q2 *cmat.FloatMatrix) {
	var qi, delta cmat.FloatMatrix

	for k := 0; k < z.Len(); k++ {
		qi.Column(Q, k)
		delta.Row(Q2, k)
		trdsecEigenVecDelta(&qi, &delta, z)
	}
}
Пример #15
0
func Sum(X *cmat.FloatMatrix, confs ...*gomas.Config) float64 {
	if X.Len() == 0 {
		return 0.0
	}
	xr, xc := X.Size()
	if xr != 1 && xc != 1 {
		return 0.0
	}
	return sum(X, X.Len())
}
Пример #16
0
func IAmax(X *cmat.FloatMatrix, confs ...*gomas.Config) int {
	if X.Len() == 0 {
		return -1
	}
	xr, xc := X.Size()
	if xr != 1 && xc != 1 {
		return -1
	}
	return iamax(X, X.Len())
}
Пример #17
0
/*
 * Triangular matrix multiplication.
 */
func MultTrm(B, A *cmat.FloatMatrix, alpha float64, bits int, confs ...*gomas.Config) *gomas.Error {
	conf := gomas.DefaultConf()
	if len(confs) > 0 {
		conf = confs[0]
	}

	if B.Len() == 0 || A.Len() == 0 {
		return nil
	}

	ok := true
	ar, ac := A.Size()
	br, bc := B.Size()
	P := ac
	E := bc
	switch {
	case bits&gomas.RIGHT != 0:
		ok = bc == ar && ar == ac
		E = br
	case bits&gomas.LEFT != 0:
		fallthrough
	default:
		ok = ac == br && ar == ac
	}
	if !ok {
		return gomas.NewError(gomas.ESIZE, "MultTrm")
	}

	// single threaded
	if conf.NProc == 1 || conf.WB <= 0 || E < conf.WB/2 {
		trmm(B, A, alpha, bits, P, 0, E, conf)
		return nil
	}

	// parallelized
	wait := make(chan int, 4)
	_, nN := blocking(0, E, conf.WB/2)
	nT := 0
	for j := 0; j < nN; j++ {
		jS := blockIndex(j, nN, conf.WB/2, E)
		jL := blockIndex(j+1, nN, conf.WB/2, E)
		task := func(q chan int) {
			trmm(B, A, alpha, bits, P, jS, jL, conf)
			q <- 1
		}
		conf.Sched.Schedule(gomas.NewTask(task, wait))
		nT += 1
	}
	for nT > 0 {
		<-wait
		nT -= 1
	}
	return nil
}
Пример #18
0
// Compute eigenvector corresponding precomputed deltas
func trdsecEigenVecDelta(qi, delta, z *cmat.FloatMatrix) {
	var dk, zk float64

	for k := 0; k < delta.Len(); k++ {
		zk = z.GetAt(k)
		dk = delta.GetAt(k)
		qi.SetAt(k, zk/dk)
	}
	s := blasd.Nrm2(qi)
	blasd.InvScale(qi, s)
}
Пример #19
0
func Nrm2(X *cmat.FloatMatrix, confs ...*gomas.Config) float64 {
	if X.Len() == 0 {
		return 0.0
	}
	xr, xc := X.Size()
	if xr != 1 && xc != 1 {
		return 0.0
	}
	nx := X.Len()
	return vnrm2(X, nx)
}
Пример #20
0
// Compute: C*(I - Y*T*Y.T) = C - C*Y*T*Y.T = C - V*T*Y.T
func updateHessRightWY(C, Y1, Y2, T, W *cmat.FloatMatrix, conf *gomas.Config) {
	var C1, C2 cmat.FloatMatrix

	if C.Len() == 0 {
		return
	}
	C1.SubMatrix(C, 0, 1, m(C), n(Y1))
	C2.SubMatrix(C, 0, 1+n(Y1), m(C), m(Y2))

	updateWithQTRight(&C1, &C2, Y1, Y2, T, W, false, conf)
}
Пример #21
0
// Compute: (I - Y*T*Y.T).T*C = C - Y*(C.T*Y*T)
func updateHessLeftWY(C, Y1, Y2, T, V *cmat.FloatMatrix, conf *gomas.Config) {
	var C1, C2 cmat.FloatMatrix

	if C.Len() == 0 {
		return
	}
	C1.SubMatrix(C, 1, 0, m(Y1), n(C))
	C2.SubMatrix(C, 1+n(Y1), 0, m(Y2), n(C))

	updateWithQTLeft(&C1, &C2, Y1, Y2, T, V, true, conf)
}
Пример #22
0
func Scale(A *cmat.FloatMatrix, alpha float64, confs ...*gomas.Config) *gomas.Error {
	if A.Len() == 0 {
		return nil
	}
	ar, ac := A.Size()
	if ar != 1 && ac != 1 {
		mscale(A, alpha, ar, ac)
		return nil
	}
	vscal(A, alpha, A.Len())
	return nil
}
Пример #23
0
/*
 * Unblocked QR decomposition with block reflector T.
 */
func unblockedQRT(A, T, W *cmat.FloatMatrix) *gomas.Error {
	var err *gomas.Error = nil
	var ATL, ATR, ABL, ABR cmat.FloatMatrix
	var A00, a10, a11, a12, A20, a21, A22 cmat.FloatMatrix
	var TTL, TTR, TBL, TBR cmat.FloatMatrix
	var T00, t01, T02, t11, t12, T22, w12 cmat.FloatMatrix

	util.Partition2x2(
		&ATL, &ATR,
		&ABL, &ABR, A, 0, 0, util.PTOPLEFT)
	util.Partition2x2(
		&TTL, &TTR,
		&TBL, &TBR, T, 0, 0, util.PTOPLEFT)

	for m(&ABR) > 0 && n(&ABR) > 0 {
		util.Repartition2x2to3x3(&ATL,
			&A00, nil, nil,
			&a10, &a11, &a12,
			&A20, &a21, &A22, A, 1, util.PBOTTOMRIGHT)
		util.Repartition2x2to3x3(&TTL,
			&T00, &t01, &T02,
			nil, &t11, &t12,
			nil, nil, &T22, T, 1, util.PBOTTOMRIGHT)

		// ------------------------------------------------------

		computeHouseholder(&a11, &a21, &t11)

		// H*[a12 A22].T
		w12.SubMatrix(W, 0, 0, a12.Len(), 1)
		applyHouseholder2x1(&t11, &a21, &a12, &A22, &w12, gomas.LEFT)

		// update T
		tauval := t11.Get(0, 0)
		if tauval != 0.0 {
			// t01 := -tauval*(a10.T + &A20.T*a21)
			//a10.CopyTo(&t01)
			blasd.Axpby(&t01, &a10, 1.0, 0.0)
			blasd.MVMult(&t01, &A20, &a21, -tauval, -tauval, gomas.TRANSA)
			// t01 := T00*t01
			blasd.MVMultTrm(&t01, &T00, 1.0, gomas.UPPER)
		}

		// ------------------------------------------------------
		util.Continue3x3to2x2(
			&ATL, &ATR,
			&ABL, &ABR, &A00, &a11, &A22, A, util.PBOTTOMRIGHT)
		util.Continue3x3to2x2(
			&TTL, &TTR,
			&TBL, &TBR, &T00, &t11, &T22, T, util.PBOTTOMRIGHT)
	}
	return err
}
Пример #24
0
func computeHouseholderVec(x, tau *cmat.FloatMatrix) {
	var alpha, x2 cmat.FloatMatrix

	r, _ := x.Size()
	alpha.SubMatrix(x, 0, 0, 1, 1)
	if r == 1 {
		x2.SubMatrix(x, 0, 1, 1, x.Len()-1) // row vector
	} else {
		x2.SubMatrix(x, 1, 0, x.Len()-1, 1)
	}
	computeHouseholder(&alpha, &x2, tau)
}
Пример #25
0
func Amax(X *cmat.FloatMatrix, confs ...*gomas.Config) float64 {
	switch X.Len() {
	case 0:
		return 0.0
	case 1:
		return X.Get(0, 0)
	}
	ix := IAmax(X, confs...)
	if ix == -1 {
		return 0.0
	}
	return X.Get(ix, 0)
}
Пример #26
0
/*
 * Compute QL factorization of a M-by-N matrix A: A = Q * L.
 *
 * Arguments:
 *  A    On entry, the M-by-N matrix A, M >= N. On exit, lower triangular matrix L
 *       and the orthogonal matrix Q as product of elementary reflectors.
 *
 * tau  On exit, the scalar factors of the elemenentary reflectors.
 *
 * W    Workspace, N-by-nb matrix used for work space in blocked invocations.
 *
 * conf The blocking configuration. If nil then default blocking configuration
 *      is used. Member conf.LB defines blocking size of blocked algorithms.
 *      If it is zero then unblocked algorithm is used.
 *
 * Returns:
 *      Error indicator.
 *
 * Additional information
 *
 *  Ortogonal matrix Q is product of elementary reflectors H(k)
 *
 *    Q = H(K-1)...H(1)H(0), where K = min(M,N)
 *
 *  Elementary reflector H(k) is stored on column k of A above the diagonal with
 *  implicit unit value on diagonal entry. The vector TAU holds scalar factors
 *  of the elementary reflectors.
 *
 *  Contents of matrix A after factorization is as follow:
 *
 *    ( v0 v1 v2 v3 )   for M=6, N=4
 *    ( v0 v1 v2 v3 )
 *    ( l  v1 v2 v3 )
 *    ( l  l  v2 v3 )
 *    ( l  l  l  v3 )
 *    ( l  l  l  l  )
 *
 *  where l is element of L, vk is element of H(k).
 *
 * DecomposeQL is compatible with lapack.DGEQLF
 */
func QLFactor(A, tau, W *cmat.FloatMatrix, confs ...*gomas.Config) *gomas.Error {
	var err *gomas.Error = nil
	var tauh cmat.FloatMatrix
	conf := gomas.CurrentConf(confs...)

	if m(A) < n(A) {
		return gomas.NewError(gomas.ESIZE, "QLFactor")
	}
	wsmin := wsQL(A, 0)
	if W == nil || W.Len() < wsmin {
		return gomas.NewError(gomas.EWORK, "QLFactor", wsmin)
	}
	if tau.Len() < n(A) {
		return gomas.NewError(gomas.ESIZE, "QLFactor")
	}
	tauh.SubMatrix(tau, 0, 0, n(A), 1)
	lb := estimateLB(A, W.Len(), wsQL)
	lb = imin(lb, conf.LB)

	if lb == 0 || n(A) <= lb {
		unblockedQL(A, &tauh, W)
	} else {
		var Twork, Wrk cmat.FloatMatrix
		// block reflector T in first LB*LB elements in workspace
		// the rest, n(A)-LB*LB, is workspace for intermediate matrix operands
		Twork.SetBuf(conf.LB, conf.LB, -1, W.Data())
		Wrk.SetBuf(n(A)-conf.LB, conf.LB, -1, W.Data()[Twork.Len():])
		blockedQL(A, &tauh, &Twork, &Wrk, lb, conf)
	}
	return err
}
Пример #27
0
// Computes eigenvectors corresponding the updated eigenvalues and rank-one update vector.
// The matrix Qd holds precomputed deltas as returned by TRDSecularSolveAll(). If Qd is nil or
// Qd same as the matrix Q then computation is in-place and Q is assumed to hold precomputed
// deltas. On exit, Q holds the column eigenvectors.
func TRDSecularEigen(Q, v, Qd *cmat.FloatMatrix, confs ...*gomas.Config) *gomas.Error {
	if m(Q) != n(Q) || (Qd != nil && (m(Qd) != n(Qd) || m(Qd) != m(Q))) {
		return gomas.NewError(gomas.ESIZE, "TRDSecularEigen")
	}
	if m(Q) != v.Len() {
		return gomas.NewError(gomas.ESIZE, "TRDSecularEigen")
	}
	if Qd == nil || Qd == Q {
		trdsecEigenBuildInplace(Q, v)
	} else {
		trdsecEigenBuild(Q, v, Qd)
	}
	return nil
}
Пример #28
0
/*
 * Unblocked code for generating M by N matrix Q with orthogonal columns which
 * are defined as the first N columns of the product of K first elementary
 * reflectors.
 *
 * Parameters nk = n(A)-K, mk = m(A)-K define the initial partitioning of
 * matrix A.
 *
 *  Q = H(k)H(k-1)...H(1)  , 0 < k <= M, where H(i) = I - tau*v*v.T
 *
 * Computation is ordered as H(k)*H(k-1)...*H(1)*I ie. from bottom to top.
 *
 * If k < M rows k+1:M are cleared and diagonal entries [k+1:M,k+1:M] are
 * set to unit. Then the matrix Q is generated by right multiplying elements below
 * of i'th elementary reflector H(i).
 *
 * Compatible to lapack.xORG2L subroutine.
 */
func unblkBuildLQ(A, Tvec, W *cmat.FloatMatrix, mk, nk int, mayClear bool) {
	var ATL, ATR, ABL, ABR cmat.FloatMatrix
	var A00, a10, a11, a12, a21, A22 cmat.FloatMatrix
	var tT, tB cmat.FloatMatrix
	var t0, tau1, t2, w12, D cmat.FloatMatrix

	util.Partition2x2(
		&ATL, &ATR,
		&ABL, &ABR, A, mk, nk, util.PBOTTOMRIGHT)
	util.Partition2x1(
		&tT,
		&tB, Tvec, mk, util.PBOTTOM)

	// zero the bottom part
	if mk > 0 && mayClear {
		blasd.Scale(&ABL, 0.0)
		blasd.Scale(&ABR, 0.0)
		D.Diag(&ABR)
		blasd.Add(&D, 1.0)
	}

	for m(&ATL) > 0 && n(&ATL) > 0 {
		util.Repartition2x2to3x3(&ATL,
			&A00, nil, nil,
			&a10, &a11, &a12,
			nil, &a21, &A22, A, 1, util.PTOPLEFT)
		util.Repartition2x1to3x1(&tT,
			&t0,
			&tau1,
			&t2, Tvec, 1, util.PTOP)
		// ------------------------------------------------------

		w12.SubMatrix(W, 0, 0, a21.Len(), 1)
		applyHouseholder2x1(&tau1, &a12, &a21, &A22, &w12, gomas.RIGHT)

		blasd.Scale(&a12, -tau1.Get(0, 0))
		a11.Set(0, 0, 1.0-tau1.Get(0, 0))

		// zero
		blasd.Scale(&a10, 0.0)
		// ------------------------------------------------------
		util.Continue3x3to2x2(
			&ATL, &ATR,
			&ABL, &ABR, &A00, &a11, &A22, A, util.PTOPLEFT)
		util.Continue3x1to2x1(
			&tT,
			&tB, &t0, &tau1, Tvec, util.PTOP)
	}
}
Пример #29
0
/*
 * Unblocked code for generating M by N matrix Q with orthogonal columns which
 * are defined as the last N columns of the product of K first elementary
 * reflectors.
 *
 * Parameter nk is last nk elementary reflectors that are not used in computing
 * the matrix Q. Parameter mk length of the first unused elementary reflectors
 * First nk columns are zeroed and subdiagonal mk-nk is set to unit.
 *
 * Compatible with lapack.DORG2L subroutine.
 */
func unblkBuildQL(A, Tvec, W *cmat.FloatMatrix, mk, nk int, mayClear bool) {
	var ATL, ATR, ABL, ABR cmat.FloatMatrix
	var A00, a01, a10, a11, a21, A22 cmat.FloatMatrix
	var tT, tB cmat.FloatMatrix
	var t0, tau1, t2, w12, D cmat.FloatMatrix

	// (mk, nk) = (rows, columns) of upper left partition
	util.Partition2x2(
		&ATL, &ATR,
		&ABL, &ABR, A, mk, nk, util.PTOPLEFT)
	util.Partition2x1(
		&tT,
		&tB, Tvec, nk, util.PTOP)

	// zero the left side
	if nk > 0 && mayClear {
		blasd.Scale(&ABL, 0.0)
		blasd.Scale(&ATL, 0.0)
		D.Diag(&ATL, nk-mk)
		blasd.Add(&D, 1.0)
	}

	for m(&ABR) > 0 && n(&ABR) > 0 {
		util.Repartition2x2to3x3(&ATL,
			&A00, &a01, nil,
			&a10, &a11, nil,
			nil, &a21, &A22, A, 1, util.PBOTTOMRIGHT)
		util.Repartition2x1to3x1(&tT,
			&t0,
			&tau1,
			&t2, Tvec, 1, util.PBOTTOM)
		// ------------------------------------------------------
		w12.SubMatrix(W, 0, 0, a10.Len(), 1)
		applyHouseholder2x1(&tau1, &a01, &a10, &A00, &w12, gomas.LEFT)

		blasd.Scale(&a01, -tau1.Get(0, 0))
		a11.Set(0, 0, 1.0-tau1.Get(0, 0))

		// zero bottom elements
		blasd.Scale(&a21, 0.0)
		// ------------------------------------------------------
		util.Continue3x3to2x2(
			&ATL, &ATR,
			&ABL, &ABR, &A00, &a11, &A22, A, util.PBOTTOMRIGHT)
		util.Continue3x1to2x1(
			&tT,
			&tB, &t0, &tau1, Tvec, util.PBOTTOM)
	}
}
Пример #30
0
func EigenSym(D, A, W *cmat.FloatMatrix, bits int, confs ...*gomas.Config) (err *gomas.Error) {

	var sD, sE, E, tau, Wred cmat.FloatMatrix
	var vv *cmat.FloatMatrix

	err = nil
	vv = nil
	conf := gomas.CurrentConf(confs...)

	if m(A) != n(A) || D.Len() != m(A) {
		err = gomas.NewError(gomas.ESIZE, "EigenSym")
		return
	}
	if bits&gomas.WANTV != 0 && W.Len() < 3*n(A) {
		err = gomas.NewError(gomas.EWORK, "EigenSym")
		return
	}

	if bits&(gomas.LOWER|gomas.UPPER) == 0 {
		bits = bits | gomas.LOWER
	}
	ioff := 1
	if bits&gomas.LOWER != 0 {
		ioff = -1
	}
	E.SetBuf(n(A)-1, 1, n(A)-1, W.Data())
	tau.SetBuf(n(A), 1, n(A), W.Data()[n(A)-1:])
	wrl := W.Len() - 2*n(A) - 1
	Wred.SetBuf(wrl, 1, wrl, W.Data()[2*n(A)-1:])

	// reduce to tridiagonal
	if err = TRDReduce(A, &tau, &Wred, bits, conf); err != nil {
		err.Update("EigenSym")
		return
	}
	sD.Diag(A)
	sE.Diag(A, ioff)
	blasd.Copy(D, &sD)
	blasd.Copy(&E, &sE)

	if bits&gomas.WANTV != 0 {
		if err = TRDBuild(A, &tau, &Wred, n(A), bits, conf); err != nil {
			err.Update("EigenSym")
			return
		}
		vv = A
	}

	// resize workspace
	wrl = W.Len() - n(A) - 1
	Wred.SetBuf(wrl, 1, wrl, W.Data()[n(A)-1:])

	if err = TRDEigen(D, &E, vv, &Wred, bits, conf); err != nil {
		err.Update("EigenSym")
		return
	}
	return
}