gonum.org/v1/gonum@v0.14.0/blas/blas64/blas64.go (about)

     1  // Copyright ©2015 The Gonum Authors. All rights reserved.
     2  // Use of this source code is governed by a BSD-style
     3  // license that can be found in the LICENSE file.
     4  
     5  package blas64
     6  
     7  import (
     8  	"gonum.org/v1/gonum/blas"
     9  	"gonum.org/v1/gonum/blas/gonum"
    10  )
    11  
    12  var blas64 blas.Float64 = gonum.Implementation{}
    13  
    14  // Use sets the BLAS float64 implementation to be used by subsequent BLAS calls.
    15  // The default implementation is
    16  // gonum.org/v1/gonum/blas/gonum.Implementation.
    17  func Use(b blas.Float64) {
    18  	blas64 = b
    19  }
    20  
    21  // Implementation returns the current BLAS float64 implementation.
    22  //
    23  // Implementation allows direct calls to the current the BLAS float64 implementation
    24  // giving finer control of parameters.
    25  func Implementation() blas.Float64 {
    26  	return blas64
    27  }
    28  
    29  // Vector represents a vector with an associated element increment.
    30  type Vector struct {
    31  	N    int
    32  	Data []float64
    33  	Inc  int
    34  }
    35  
    36  // General represents a matrix using the conventional storage scheme.
    37  type General struct {
    38  	Rows, Cols int
    39  	Data       []float64
    40  	Stride     int
    41  }
    42  
    43  // Band represents a band matrix using the band storage scheme.
    44  type Band struct {
    45  	Rows, Cols int
    46  	KL, KU     int
    47  	Data       []float64
    48  	Stride     int
    49  }
    50  
    51  // Triangular represents a triangular matrix using the conventional storage scheme.
    52  type Triangular struct {
    53  	Uplo   blas.Uplo
    54  	Diag   blas.Diag
    55  	N      int
    56  	Data   []float64
    57  	Stride int
    58  }
    59  
    60  // TriangularBand represents a triangular matrix using the band storage scheme.
    61  type TriangularBand struct {
    62  	Uplo   blas.Uplo
    63  	Diag   blas.Diag
    64  	N, K   int
    65  	Data   []float64
    66  	Stride int
    67  }
    68  
    69  // TriangularPacked represents a triangular matrix using the packed storage scheme.
    70  type TriangularPacked struct {
    71  	Uplo blas.Uplo
    72  	Diag blas.Diag
    73  	N    int
    74  	Data []float64
    75  }
    76  
    77  // Symmetric represents a symmetric matrix using the conventional storage scheme.
    78  type Symmetric struct {
    79  	Uplo   blas.Uplo
    80  	N      int
    81  	Data   []float64
    82  	Stride int
    83  }
    84  
    85  // SymmetricBand represents a symmetric matrix using the band storage scheme.
    86  type SymmetricBand struct {
    87  	Uplo   blas.Uplo
    88  	N, K   int
    89  	Data   []float64
    90  	Stride int
    91  }
    92  
    93  // SymmetricPacked represents a symmetric matrix using the packed storage scheme.
    94  type SymmetricPacked struct {
    95  	Uplo blas.Uplo
    96  	N    int
    97  	Data []float64
    98  }
    99  
   100  // Level 1
   101  
   102  const (
   103  	negInc    = "blas64: negative vector increment"
   104  	badLength = "blas64: vector length mismatch"
   105  )
   106  
   107  // Dot computes the dot product of the two vectors:
   108  //
   109  //	\sum_i x[i]*y[i].
   110  //
   111  // Dot will panic if the lengths of x and y do not match.
   112  func Dot(x, y Vector) float64 {
   113  	if x.N != y.N {
   114  		panic(badLength)
   115  	}
   116  	return blas64.Ddot(x.N, x.Data, x.Inc, y.Data, y.Inc)
   117  }
   118  
   119  // Nrm2 computes the Euclidean norm of the vector x:
   120  //
   121  //	sqrt(\sum_i x[i]*x[i]).
   122  //
   123  // Nrm2 will panic if the vector increment is negative.
   124  func Nrm2(x Vector) float64 {
   125  	if x.Inc < 0 {
   126  		panic(negInc)
   127  	}
   128  	return blas64.Dnrm2(x.N, x.Data, x.Inc)
   129  }
   130  
   131  // Asum computes the sum of the absolute values of the elements of x:
   132  //
   133  //	\sum_i |x[i]|.
   134  //
   135  // Asum will panic if the vector increment is negative.
   136  func Asum(x Vector) float64 {
   137  	if x.Inc < 0 {
   138  		panic(negInc)
   139  	}
   140  	return blas64.Dasum(x.N, x.Data, x.Inc)
   141  }
   142  
   143  // Iamax returns the index of an element of x with the largest absolute value.
   144  // If there are multiple such indices the earliest is returned.
   145  // Iamax returns -1 if n == 0.
   146  //
   147  // Iamax will panic if the vector increment is negative.
   148  func Iamax(x Vector) int {
   149  	if x.Inc < 0 {
   150  		panic(negInc)
   151  	}
   152  	return blas64.Idamax(x.N, x.Data, x.Inc)
   153  }
   154  
   155  // Swap exchanges the elements of the two vectors:
   156  //
   157  //	x[i], y[i] = y[i], x[i] for all i.
   158  //
   159  // Swap will panic if the lengths of x and y do not match.
   160  func Swap(x, y Vector) {
   161  	if x.N != y.N {
   162  		panic(badLength)
   163  	}
   164  	blas64.Dswap(x.N, x.Data, x.Inc, y.Data, y.Inc)
   165  }
   166  
   167  // Copy copies the elements of x into the elements of y:
   168  //
   169  //	y[i] = x[i] for all i.
   170  //
   171  // Copy will panic if the lengths of x and y do not match.
   172  func Copy(x, y Vector) {
   173  	if x.N != y.N {
   174  		panic(badLength)
   175  	}
   176  	blas64.Dcopy(x.N, x.Data, x.Inc, y.Data, y.Inc)
   177  }
   178  
   179  // Axpy adds x scaled by alpha to y:
   180  //
   181  //	y[i] += alpha*x[i] for all i.
   182  //
   183  // Axpy will panic if the lengths of x and y do not match.
   184  func Axpy(alpha float64, x, y Vector) {
   185  	if x.N != y.N {
   186  		panic(badLength)
   187  	}
   188  	blas64.Daxpy(x.N, alpha, x.Data, x.Inc, y.Data, y.Inc)
   189  }
   190  
   191  // Rotg computes the parameters of a Givens plane rotation so that
   192  //
   193  //	⎡ c s⎤   ⎡a⎤   ⎡r⎤
   194  //	⎣-s c⎦ * ⎣b⎦ = ⎣0⎦
   195  //
   196  // where a and b are the Cartesian coordinates of a given point.
   197  // c, s, and r are defined as
   198  //
   199  //	r = ±Sqrt(a^2 + b^2),
   200  //	c = a/r, the cosine of the rotation angle,
   201  //	s = a/r, the sine of the rotation angle,
   202  //
   203  // and z is defined such that
   204  //
   205  //	if |a| > |b|,        z = s,
   206  //	otherwise if c != 0, z = 1/c,
   207  //	otherwise            z = 1.
   208  func Rotg(a, b float64) (c, s, r, z float64) {
   209  	return blas64.Drotg(a, b)
   210  }
   211  
   212  // Rotmg computes the modified Givens rotation. See
   213  // http://www.netlib.org/lapack/explore-html/df/deb/drotmg_8f.html
   214  // for more details.
   215  func Rotmg(d1, d2, b1, b2 float64) (p blas.DrotmParams, rd1, rd2, rb1 float64) {
   216  	return blas64.Drotmg(d1, d2, b1, b2)
   217  }
   218  
   219  // Rot applies a plane transformation to n points represented by the vectors x
   220  // and y:
   221  //
   222  //	x[i] =  c*x[i] + s*y[i],
   223  //	y[i] = -s*x[i] + c*y[i], for all i.
   224  func Rot(x, y Vector, c, s float64) {
   225  	if x.N != y.N {
   226  		panic(badLength)
   227  	}
   228  	blas64.Drot(x.N, x.Data, x.Inc, y.Data, y.Inc, c, s)
   229  }
   230  
   231  // Rotm applies the modified Givens rotation to n points represented by the
   232  // vectors x and y.
   233  func Rotm(x, y Vector, p blas.DrotmParams) {
   234  	if x.N != y.N {
   235  		panic(badLength)
   236  	}
   237  	blas64.Drotm(x.N, x.Data, x.Inc, y.Data, y.Inc, p)
   238  }
   239  
   240  // Scal scales the vector x by alpha:
   241  //
   242  //	x[i] *= alpha for all i.
   243  //
   244  // Scal will panic if the vector increment is negative.
   245  func Scal(alpha float64, x Vector) {
   246  	if x.Inc < 0 {
   247  		panic(negInc)
   248  	}
   249  	blas64.Dscal(x.N, alpha, x.Data, x.Inc)
   250  }
   251  
   252  // Level 2
   253  
   254  // Gemv computes
   255  //
   256  //	y = alpha * A * x + beta * y   if t == blas.NoTrans,
   257  //	y = alpha * Aᵀ * x + beta * y  if t == blas.Trans or blas.ConjTrans,
   258  //
   259  // where A is an m×n dense matrix, x and y are vectors, and alpha and beta are scalars.
   260  func Gemv(t blas.Transpose, alpha float64, a General, x Vector, beta float64, y Vector) {
   261  	blas64.Dgemv(t, a.Rows, a.Cols, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   262  }
   263  
   264  // Gbmv computes
   265  //
   266  //	y = alpha * A * x + beta * y   if t == blas.NoTrans,
   267  //	y = alpha * Aᵀ * x + beta * y  if t == blas.Trans or blas.ConjTrans,
   268  //
   269  // where A is an m×n band matrix, x and y are vectors, and alpha and beta are scalars.
   270  func Gbmv(t blas.Transpose, alpha float64, a Band, x Vector, beta float64, y Vector) {
   271  	blas64.Dgbmv(t, a.Rows, a.Cols, a.KL, a.KU, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   272  }
   273  
   274  // Trmv computes
   275  //
   276  //	x = A * x   if t == blas.NoTrans,
   277  //	x = Aᵀ * x  if t == blas.Trans or blas.ConjTrans,
   278  //
   279  // where A is an n×n triangular matrix, and x is a vector.
   280  func Trmv(t blas.Transpose, a Triangular, x Vector) {
   281  	blas64.Dtrmv(a.Uplo, t, a.Diag, a.N, a.Data, a.Stride, x.Data, x.Inc)
   282  }
   283  
   284  // Tbmv computes
   285  //
   286  //	x = A * x   if t == blas.NoTrans,
   287  //	x = Aᵀ * x  if t == blas.Trans or blas.ConjTrans,
   288  //
   289  // where A is an n×n triangular band matrix, and x is a vector.
   290  func Tbmv(t blas.Transpose, a TriangularBand, x Vector) {
   291  	blas64.Dtbmv(a.Uplo, t, a.Diag, a.N, a.K, a.Data, a.Stride, x.Data, x.Inc)
   292  }
   293  
   294  // Tpmv computes
   295  //
   296  //	x = A * x   if t == blas.NoTrans,
   297  //	x = Aᵀ * x  if t == blas.Trans or blas.ConjTrans,
   298  //
   299  // where A is an n×n triangular matrix in packed format, and x is a vector.
   300  func Tpmv(t blas.Transpose, a TriangularPacked, x Vector) {
   301  	blas64.Dtpmv(a.Uplo, t, a.Diag, a.N, a.Data, x.Data, x.Inc)
   302  }
   303  
   304  // Trsv solves
   305  //
   306  //	A * x = b   if t == blas.NoTrans,
   307  //	Aᵀ * x = b  if t == blas.Trans or blas.ConjTrans,
   308  //
   309  // where A is an n×n triangular matrix, and x and b are vectors.
   310  //
   311  // At entry to the function, x contains the values of b, and the result is
   312  // stored in-place into x.
   313  //
   314  // No test for singularity or near-singularity is included in this
   315  // routine. Such tests must be performed before calling this routine.
   316  func Trsv(t blas.Transpose, a Triangular, x Vector) {
   317  	blas64.Dtrsv(a.Uplo, t, a.Diag, a.N, a.Data, a.Stride, x.Data, x.Inc)
   318  }
   319  
   320  // Tbsv solves
   321  //
   322  //	A * x = b   if t == blas.NoTrans,
   323  //	Aᵀ * x = b  if t == blas.Trans or blas.ConjTrans,
   324  //
   325  // where A is an n×n triangular band matrix, and x and b are vectors.
   326  //
   327  // At entry to the function, x contains the values of b, and the result is
   328  // stored in place into x.
   329  //
   330  // No test for singularity or near-singularity is included in this
   331  // routine. Such tests must be performed before calling this routine.
   332  func Tbsv(t blas.Transpose, a TriangularBand, x Vector) {
   333  	blas64.Dtbsv(a.Uplo, t, a.Diag, a.N, a.K, a.Data, a.Stride, x.Data, x.Inc)
   334  }
   335  
   336  // Tpsv solves
   337  //
   338  //	A * x = b   if t == blas.NoTrans,
   339  //	Aᵀ * x = b  if t == blas.Trans or blas.ConjTrans,
   340  //
   341  // where A is an n×n triangular matrix in packed format, and x and b are
   342  // vectors.
   343  //
   344  // At entry to the function, x contains the values of b, and the result is
   345  // stored in place into x.
   346  //
   347  // No test for singularity or near-singularity is included in this
   348  // routine. Such tests must be performed before calling this routine.
   349  func Tpsv(t blas.Transpose, a TriangularPacked, x Vector) {
   350  	blas64.Dtpsv(a.Uplo, t, a.Diag, a.N, a.Data, x.Data, x.Inc)
   351  }
   352  
   353  // Symv computes
   354  //
   355  //	y = alpha * A * x + beta * y,
   356  //
   357  // where A is an n×n symmetric matrix, x and y are vectors, and alpha and
   358  // beta are scalars.
   359  func Symv(alpha float64, a Symmetric, x Vector, beta float64, y Vector) {
   360  	blas64.Dsymv(a.Uplo, a.N, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   361  }
   362  
   363  // Sbmv performs
   364  //
   365  //	y = alpha * A * x + beta * y,
   366  //
   367  // where A is an n×n symmetric band matrix, x and y are vectors, and alpha
   368  // and beta are scalars.
   369  func Sbmv(alpha float64, a SymmetricBand, x Vector, beta float64, y Vector) {
   370  	blas64.Dsbmv(a.Uplo, a.N, a.K, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   371  }
   372  
   373  // Spmv performs
   374  //
   375  //	y = alpha * A * x + beta * y,
   376  //
   377  // where A is an n×n symmetric matrix in packed format, x and y are vectors,
   378  // and alpha and beta are scalars.
   379  func Spmv(alpha float64, a SymmetricPacked, x Vector, beta float64, y Vector) {
   380  	blas64.Dspmv(a.Uplo, a.N, alpha, a.Data, x.Data, x.Inc, beta, y.Data, y.Inc)
   381  }
   382  
   383  // Ger performs a rank-1 update
   384  //
   385  //	A += alpha * x * yᵀ,
   386  //
   387  // where A is an m×n dense matrix, x and y are vectors, and alpha is a scalar.
   388  func Ger(alpha float64, x, y Vector, a General) {
   389  	blas64.Dger(a.Rows, a.Cols, alpha, x.Data, x.Inc, y.Data, y.Inc, a.Data, a.Stride)
   390  }
   391  
   392  // Syr performs a rank-1 update
   393  //
   394  //	A += alpha * x * xᵀ,
   395  //
   396  // where A is an n×n symmetric matrix, x is a vector, and alpha is a scalar.
   397  func Syr(alpha float64, x Vector, a Symmetric) {
   398  	blas64.Dsyr(a.Uplo, a.N, alpha, x.Data, x.Inc, a.Data, a.Stride)
   399  }
   400  
   401  // Spr performs the rank-1 update
   402  //
   403  //	A += alpha * x * xᵀ,
   404  //
   405  // where A is an n×n symmetric matrix in packed format, x is a vector, and
   406  // alpha is a scalar.
   407  func Spr(alpha float64, x Vector, a SymmetricPacked) {
   408  	blas64.Dspr(a.Uplo, a.N, alpha, x.Data, x.Inc, a.Data)
   409  }
   410  
   411  // Syr2 performs a rank-2 update
   412  //
   413  //	A += alpha * x * yᵀ + alpha * y * xᵀ,
   414  //
   415  // where A is a symmetric n×n matrix, x and y are vectors, and alpha is a scalar.
   416  func Syr2(alpha float64, x, y Vector, a Symmetric) {
   417  	blas64.Dsyr2(a.Uplo, a.N, alpha, x.Data, x.Inc, y.Data, y.Inc, a.Data, a.Stride)
   418  }
   419  
   420  // Spr2 performs a rank-2 update
   421  //
   422  //	A += alpha * x * yᵀ + alpha * y * xᵀ,
   423  //
   424  // where A is an n×n symmetric matrix in packed format, x and y are vectors,
   425  // and alpha is a scalar.
   426  func Spr2(alpha float64, x, y Vector, a SymmetricPacked) {
   427  	blas64.Dspr2(a.Uplo, a.N, alpha, x.Data, x.Inc, y.Data, y.Inc, a.Data)
   428  }
   429  
   430  // Level 3
   431  
   432  // Gemm computes
   433  //
   434  //	C = alpha * A * B + beta * C,
   435  //
   436  // where A, B, and C are dense matrices, and alpha and beta are scalars.
   437  // tA and tB specify whether A or B are transposed.
   438  func Gemm(tA, tB blas.Transpose, alpha float64, a, b General, beta float64, c General) {
   439  	var m, n, k int
   440  	if tA == blas.NoTrans {
   441  		m, k = a.Rows, a.Cols
   442  	} else {
   443  		m, k = a.Cols, a.Rows
   444  	}
   445  	if tB == blas.NoTrans {
   446  		n = b.Cols
   447  	} else {
   448  		n = b.Rows
   449  	}
   450  	blas64.Dgemm(tA, tB, m, n, k, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   451  }
   452  
   453  // Symm performs
   454  //
   455  //	C = alpha * A * B + beta * C  if s == blas.Left,
   456  //	C = alpha * B * A + beta * C  if s == blas.Right,
   457  //
   458  // where A is an n×n or m×m symmetric matrix, B and C are m×n matrices, and
   459  // alpha is a scalar.
   460  func Symm(s blas.Side, alpha float64, a Symmetric, b General, beta float64, c General) {
   461  	var m, n int
   462  	if s == blas.Left {
   463  		m, n = a.N, b.Cols
   464  	} else {
   465  		m, n = b.Rows, a.N
   466  	}
   467  	blas64.Dsymm(s, a.Uplo, m, n, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   468  }
   469  
   470  // Syrk performs a symmetric rank-k update
   471  //
   472  //	C = alpha * A * Aᵀ + beta * C  if t == blas.NoTrans,
   473  //	C = alpha * Aᵀ * A + beta * C  if t == blas.Trans or blas.ConjTrans,
   474  //
   475  // where C is an n×n symmetric matrix, A is an n×k matrix if t == blas.NoTrans and
   476  // a k×n matrix otherwise, and alpha and beta are scalars.
   477  func Syrk(t blas.Transpose, alpha float64, a General, beta float64, c Symmetric) {
   478  	var n, k int
   479  	if t == blas.NoTrans {
   480  		n, k = a.Rows, a.Cols
   481  	} else {
   482  		n, k = a.Cols, a.Rows
   483  	}
   484  	blas64.Dsyrk(c.Uplo, t, n, k, alpha, a.Data, a.Stride, beta, c.Data, c.Stride)
   485  }
   486  
   487  // Syr2k performs a symmetric rank-2k update
   488  //
   489  //	C = alpha * A * Bᵀ + alpha * B * Aᵀ + beta * C  if t == blas.NoTrans,
   490  //	C = alpha * Aᵀ * B + alpha * Bᵀ * A + beta * C  if t == blas.Trans or blas.ConjTrans,
   491  //
   492  // where C is an n×n symmetric matrix, A and B are n×k matrices if t == NoTrans
   493  // and k×n matrices otherwise, and alpha and beta are scalars.
   494  func Syr2k(t blas.Transpose, alpha float64, a, b General, beta float64, c Symmetric) {
   495  	var n, k int
   496  	if t == blas.NoTrans {
   497  		n, k = a.Rows, a.Cols
   498  	} else {
   499  		n, k = a.Cols, a.Rows
   500  	}
   501  	blas64.Dsyr2k(c.Uplo, t, n, k, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   502  }
   503  
   504  // Trmm performs
   505  //
   506  //	B = alpha * A * B   if tA == blas.NoTrans and s == blas.Left,
   507  //	B = alpha * Aᵀ * B  if tA == blas.Trans or blas.ConjTrans, and s == blas.Left,
   508  //	B = alpha * B * A   if tA == blas.NoTrans and s == blas.Right,
   509  //	B = alpha * B * Aᵀ  if tA == blas.Trans or blas.ConjTrans, and s == blas.Right,
   510  //
   511  // where A is an n×n or m×m triangular matrix, B is an m×n matrix, and alpha is
   512  // a scalar.
   513  func Trmm(s blas.Side, tA blas.Transpose, alpha float64, a Triangular, b General) {
   514  	blas64.Dtrmm(s, a.Uplo, tA, a.Diag, b.Rows, b.Cols, alpha, a.Data, a.Stride, b.Data, b.Stride)
   515  }
   516  
   517  // Trsm solves
   518  //
   519  //	A * X = alpha * B   if tA == blas.NoTrans and s == blas.Left,
   520  //	Aᵀ * X = alpha * B  if tA == blas.Trans or blas.ConjTrans, and s == blas.Left,
   521  //	X * A = alpha * B   if tA == blas.NoTrans and s == blas.Right,
   522  //	X * Aᵀ = alpha * B  if tA == blas.Trans or blas.ConjTrans, and s == blas.Right,
   523  //
   524  // where A is an n×n or m×m triangular matrix, X and B are m×n matrices, and
   525  // alpha is a scalar.
   526  //
   527  // At entry to the function, X contains the values of B, and the result is
   528  // stored in-place into X.
   529  //
   530  // No check is made that A is invertible.
   531  func Trsm(s blas.Side, tA blas.Transpose, alpha float64, a Triangular, b General) {
   532  	blas64.Dtrsm(s, a.Uplo, tA, a.Diag, b.Rows, b.Cols, alpha, a.Data, a.Stride, b.Data, b.Stride)
   533  }