gonum.org/v1/gonum@v0.14.0/blas/cblas64/cblas64.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 cblas64
     6  
     7  import (
     8  	"gonum.org/v1/gonum/blas"
     9  	"gonum.org/v1/gonum/blas/gonum"
    10  )
    11  
    12  var cblas64 blas.Complex64 = gonum.Implementation{}
    13  
    14  // Use sets the BLAS complex64 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.Complex64) {
    18  	cblas64 = b
    19  }
    20  
    21  // Implementation returns the current BLAS complex64 implementation.
    22  //
    23  // Implementation allows direct calls to the current the BLAS complex64 implementation
    24  // giving finer control of parameters.
    25  func Implementation() blas.Complex64 {
    26  	return cblas64
    27  }
    28  
    29  // Vector represents a vector with an associated element increment.
    30  type Vector struct {
    31  	N    int
    32  	Inc  int
    33  	Data []complex64
    34  }
    35  
    36  // General represents a matrix using the conventional storage scheme.
    37  type General struct {
    38  	Rows, Cols int
    39  	Stride     int
    40  	Data       []complex64
    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  	Stride     int
    48  	Data       []complex64
    49  }
    50  
    51  // Triangular represents a triangular matrix using the conventional storage scheme.
    52  type Triangular struct {
    53  	N      int
    54  	Stride int
    55  	Data   []complex64
    56  	Uplo   blas.Uplo
    57  	Diag   blas.Diag
    58  }
    59  
    60  // TriangularBand represents a triangular matrix using the band storage scheme.
    61  type TriangularBand struct {
    62  	N, K   int
    63  	Stride int
    64  	Data   []complex64
    65  	Uplo   blas.Uplo
    66  	Diag   blas.Diag
    67  }
    68  
    69  // TriangularPacked represents a triangular matrix using the packed storage scheme.
    70  type TriangularPacked struct {
    71  	N    int
    72  	Data []complex64
    73  	Uplo blas.Uplo
    74  	Diag blas.Diag
    75  }
    76  
    77  // Symmetric represents a symmetric matrix using the conventional storage scheme.
    78  type Symmetric struct {
    79  	N      int
    80  	Stride int
    81  	Data   []complex64
    82  	Uplo   blas.Uplo
    83  }
    84  
    85  // SymmetricBand represents a symmetric matrix using the band storage scheme.
    86  type SymmetricBand struct {
    87  	N, K   int
    88  	Stride int
    89  	Data   []complex64
    90  	Uplo   blas.Uplo
    91  }
    92  
    93  // SymmetricPacked represents a symmetric matrix using the packed storage scheme.
    94  type SymmetricPacked struct {
    95  	N    int
    96  	Data []complex64
    97  	Uplo blas.Uplo
    98  }
    99  
   100  // Hermitian represents an Hermitian matrix using the conventional storage scheme.
   101  type Hermitian Symmetric
   102  
   103  // HermitianBand represents an Hermitian matrix using the band storage scheme.
   104  type HermitianBand SymmetricBand
   105  
   106  // HermitianPacked represents an Hermitian matrix using the packed storage scheme.
   107  type HermitianPacked SymmetricPacked
   108  
   109  // Level 1
   110  
   111  const (
   112  	negInc    = "cblas64: negative vector increment"
   113  	badLength = "cblas64: vector length mismatch"
   114  )
   115  
   116  // Dotu computes the dot product of the two vectors without
   117  // complex conjugation:
   118  //
   119  //	xᵀ * y
   120  //
   121  // Dotu will panic if the lengths of x and y do not match.
   122  func Dotu(x, y Vector) complex64 {
   123  	if x.N != y.N {
   124  		panic(badLength)
   125  	}
   126  	return cblas64.Cdotu(x.N, x.Data, x.Inc, y.Data, y.Inc)
   127  }
   128  
   129  // Dotc computes the dot product of the two vectors with
   130  // complex conjugation:
   131  //
   132  //	xᴴ * y.
   133  //
   134  // Dotc will panic if the lengths of x and y do not match.
   135  func Dotc(x, y Vector) complex64 {
   136  	if x.N != y.N {
   137  		panic(badLength)
   138  	}
   139  	return cblas64.Cdotc(x.N, x.Data, x.Inc, y.Data, y.Inc)
   140  }
   141  
   142  // Nrm2 computes the Euclidean norm of the vector x:
   143  //
   144  //	sqrt(\sum_i x[i] * x[i]).
   145  //
   146  // Nrm2 will panic if the vector increment is negative.
   147  func Nrm2(x Vector) float32 {
   148  	if x.Inc < 0 {
   149  		panic(negInc)
   150  	}
   151  	return cblas64.Scnrm2(x.N, x.Data, x.Inc)
   152  }
   153  
   154  // Asum computes the sum of magnitudes of the real and imaginary parts of
   155  // elements of the vector x:
   156  //
   157  //	\sum_i (|Re x[i]| + |Im x[i]|).
   158  //
   159  // Asum will panic if the vector increment is negative.
   160  func Asum(x Vector) float32 {
   161  	if x.Inc < 0 {
   162  		panic(negInc)
   163  	}
   164  	return cblas64.Scasum(x.N, x.Data, x.Inc)
   165  }
   166  
   167  // Iamax returns the index of an element of x with the largest sum of
   168  // magnitudes of the real and imaginary parts (|Re x[i]|+|Im x[i]|).
   169  // If there are multiple such indices, the earliest is returned.
   170  //
   171  // Iamax returns -1 if n == 0.
   172  //
   173  // Iamax will panic if the vector increment is negative.
   174  func Iamax(x Vector) int {
   175  	if x.Inc < 0 {
   176  		panic(negInc)
   177  	}
   178  	return cblas64.Icamax(x.N, x.Data, x.Inc)
   179  }
   180  
   181  // Swap exchanges the elements of two vectors:
   182  //
   183  //	x[i], y[i] = y[i], x[i] for all i.
   184  //
   185  // Swap will panic if the lengths of x and y do not match.
   186  func Swap(x, y Vector) {
   187  	if x.N != y.N {
   188  		panic(badLength)
   189  	}
   190  	cblas64.Cswap(x.N, x.Data, x.Inc, y.Data, y.Inc)
   191  }
   192  
   193  // Copy copies the elements of x into the elements of y:
   194  //
   195  //	y[i] = x[i] for all i.
   196  //
   197  // Copy will panic if the lengths of x and y do not match.
   198  func Copy(x, y Vector) {
   199  	if x.N != y.N {
   200  		panic(badLength)
   201  	}
   202  	cblas64.Ccopy(x.N, x.Data, x.Inc, y.Data, y.Inc)
   203  }
   204  
   205  // Axpy computes
   206  //
   207  //	y = alpha * x + y,
   208  //
   209  // where x and y are vectors, and alpha is a scalar.
   210  // Axpy will panic if the lengths of x and y do not match.
   211  func Axpy(alpha complex64, x, y Vector) {
   212  	if x.N != y.N {
   213  		panic(badLength)
   214  	}
   215  	cblas64.Caxpy(x.N, alpha, x.Data, x.Inc, y.Data, y.Inc)
   216  }
   217  
   218  // Scal computes
   219  //
   220  //	x = alpha * x,
   221  //
   222  // where x is a vector, and alpha is a scalar.
   223  //
   224  // Scal will panic if the vector increment is negative.
   225  func Scal(alpha complex64, x Vector) {
   226  	if x.Inc < 0 {
   227  		panic(negInc)
   228  	}
   229  	cblas64.Cscal(x.N, alpha, x.Data, x.Inc)
   230  }
   231  
   232  // Dscal computes
   233  //
   234  //	x = alpha * x,
   235  //
   236  // where x is a vector, and alpha is a real scalar.
   237  //
   238  // Dscal will panic if the vector increment is negative.
   239  func Dscal(alpha float32, x Vector) {
   240  	if x.Inc < 0 {
   241  		panic(negInc)
   242  	}
   243  	cblas64.Csscal(x.N, alpha, x.Data, x.Inc)
   244  }
   245  
   246  // Level 2
   247  
   248  // Gemv computes
   249  //
   250  //	y = alpha * A * x + beta * y   if t == blas.NoTrans,
   251  //	y = alpha * Aᵀ * x + beta * y  if t == blas.Trans,
   252  //	y = alpha * Aᴴ * x + beta * y  if t == blas.ConjTrans,
   253  //
   254  // where A is an m×n dense matrix, x and y are vectors, and alpha and beta are
   255  // scalars.
   256  func Gemv(t blas.Transpose, alpha complex64, a General, x Vector, beta complex64, y Vector) {
   257  	cblas64.Cgemv(t, a.Rows, a.Cols, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   258  }
   259  
   260  // Gbmv computes
   261  //
   262  //	y = alpha * A * x + beta * y   if t == blas.NoTrans,
   263  //	y = alpha * Aᵀ * x + beta * y  if t == blas.Trans,
   264  //	y = alpha * Aᴴ * x + beta * y  if t == blas.ConjTrans,
   265  //
   266  // where A is an m×n band matrix, x and y are vectors, and alpha and beta are
   267  // scalars.
   268  func Gbmv(t blas.Transpose, alpha complex64, a Band, x Vector, beta complex64, y Vector) {
   269  	cblas64.Cgbmv(t, a.Rows, a.Cols, a.KL, a.KU, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   270  }
   271  
   272  // Trmv computes
   273  //
   274  //	x = A * x   if t == blas.NoTrans,
   275  //	x = Aᵀ * x  if t == blas.Trans,
   276  //	x = Aᴴ * x  if t == blas.ConjTrans,
   277  //
   278  // where A is an n×n triangular matrix, and x is a vector.
   279  func Trmv(t blas.Transpose, a Triangular, x Vector) {
   280  	cblas64.Ctrmv(a.Uplo, t, a.Diag, a.N, a.Data, a.Stride, x.Data, x.Inc)
   281  }
   282  
   283  // Tbmv computes
   284  //
   285  //	x = A * x   if t == blas.NoTrans,
   286  //	x = Aᵀ * x  if t == blas.Trans,
   287  //	x = Aᴴ * x  if t == 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  	cblas64.Ctbmv(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,
   298  //	x = Aᴴ * x  if t == blas.ConjTrans,
   299  //
   300  // where A is an n×n triangular matrix in packed format, and x is a vector.
   301  func Tpmv(t blas.Transpose, a TriangularPacked, x Vector) {
   302  	cblas64.Ctpmv(a.Uplo, t, a.Diag, a.N, a.Data, x.Data, x.Inc)
   303  }
   304  
   305  // Trsv solves
   306  //
   307  //	A * x = b   if t == blas.NoTrans,
   308  //	Aᵀ * x = b  if t == blas.Trans,
   309  //	Aᴴ * x = b  if t == blas.ConjTrans,
   310  //
   311  // where A is an n×n triangular matrix and x is a vector.
   312  //
   313  // At entry to the function, x contains the values of b, and the result is
   314  // stored in-place into x.
   315  //
   316  // No test for singularity or near-singularity is included in this
   317  // routine. Such tests must be performed before calling this routine.
   318  func Trsv(t blas.Transpose, a Triangular, x Vector) {
   319  	cblas64.Ctrsv(a.Uplo, t, a.Diag, a.N, a.Data, a.Stride, x.Data, x.Inc)
   320  }
   321  
   322  // Tbsv solves
   323  //
   324  //	A * x = b   if t == blas.NoTrans,
   325  //	Aᵀ * x = b  if t == blas.Trans,
   326  //	Aᴴ * x = b  if t == blas.ConjTrans,
   327  //
   328  // where A is an n×n triangular band matrix, and x is a vector.
   329  //
   330  // At entry to the function, x contains the values of b, and the result is
   331  // stored in-place into x.
   332  //
   333  // No test for singularity or near-singularity is included in this
   334  // routine. Such tests must be performed before calling this routine.
   335  func Tbsv(t blas.Transpose, a TriangularBand, x Vector) {
   336  	cblas64.Ctbsv(a.Uplo, t, a.Diag, a.N, a.K, a.Data, a.Stride, x.Data, x.Inc)
   337  }
   338  
   339  // Tpsv solves
   340  //
   341  //	A * x = b   if t == blas.NoTrans,
   342  //	Aᵀ * x = b  if t == blas.Trans,
   343  //	Aᴴ * x = b  if t == blas.ConjTrans,
   344  //
   345  // where A is an n×n triangular matrix in packed format and x is a vector.
   346  //
   347  // At entry to the function, x contains the values of b, and the result is
   348  // stored in-place into x.
   349  //
   350  // No test for singularity or near-singularity is included in this
   351  // routine. Such tests must be performed before calling this routine.
   352  func Tpsv(t blas.Transpose, a TriangularPacked, x Vector) {
   353  	cblas64.Ctpsv(a.Uplo, t, a.Diag, a.N, a.Data, x.Data, x.Inc)
   354  }
   355  
   356  // Hemv computes
   357  //
   358  //	y = alpha * A * x + beta * y,
   359  //
   360  // where A is an n×n Hermitian matrix, x and y are vectors, and alpha and
   361  // beta are scalars.
   362  func Hemv(alpha complex64, a Hermitian, x Vector, beta complex64, y Vector) {
   363  	cblas64.Chemv(a.Uplo, a.N, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   364  }
   365  
   366  // Hbmv performs
   367  //
   368  //	y = alpha * A * x + beta * y,
   369  //
   370  // where A is an n×n Hermitian band matrix, x and y are vectors, and alpha
   371  // and beta are scalars.
   372  func Hbmv(alpha complex64, a HermitianBand, x Vector, beta complex64, y Vector) {
   373  	cblas64.Chbmv(a.Uplo, a.N, a.K, alpha, a.Data, a.Stride, x.Data, x.Inc, beta, y.Data, y.Inc)
   374  }
   375  
   376  // Hpmv performs
   377  //
   378  //	y = alpha * A * x + beta * y,
   379  //
   380  // where A is an n×n Hermitian matrix in packed format, x and y are vectors,
   381  // and alpha and beta are scalars.
   382  func Hpmv(alpha complex64, a HermitianPacked, x Vector, beta complex64, y Vector) {
   383  	cblas64.Chpmv(a.Uplo, a.N, alpha, a.Data, x.Data, x.Inc, beta, y.Data, y.Inc)
   384  }
   385  
   386  // Geru performs a rank-1 update
   387  //
   388  //	A += alpha * x * yᵀ,
   389  //
   390  // where A is an m×n dense matrix, x and y are vectors, and alpha is a scalar.
   391  func Geru(alpha complex64, x, y Vector, a General) {
   392  	cblas64.Cgeru(a.Rows, a.Cols, alpha, x.Data, x.Inc, y.Data, y.Inc, a.Data, a.Stride)
   393  }
   394  
   395  // Gerc performs a rank-1 update
   396  //
   397  //	A += alpha * x * yᴴ,
   398  //
   399  // where A is an m×n dense matrix, x and y are vectors, and alpha is a scalar.
   400  func Gerc(alpha complex64, x, y Vector, a General) {
   401  	cblas64.Cgerc(a.Rows, a.Cols, alpha, x.Data, x.Inc, y.Data, y.Inc, a.Data, a.Stride)
   402  }
   403  
   404  // Her performs a rank-1 update
   405  //
   406  //	A += alpha * x * yᵀ,
   407  //
   408  // where A is an m×n Hermitian matrix, x and y are vectors, and alpha is a scalar.
   409  func Her(alpha float32, x Vector, a Hermitian) {
   410  	cblas64.Cher(a.Uplo, a.N, alpha, x.Data, x.Inc, a.Data, a.Stride)
   411  }
   412  
   413  // Hpr performs a rank-1 update
   414  //
   415  //	A += alpha * x * xᴴ,
   416  //
   417  // where A is an n×n Hermitian matrix in packed format, x is a vector, and
   418  // alpha is a scalar.
   419  func Hpr(alpha float32, x Vector, a HermitianPacked) {
   420  	cblas64.Chpr(a.Uplo, a.N, alpha, x.Data, x.Inc, a.Data)
   421  }
   422  
   423  // Her2 performs a rank-2 update
   424  //
   425  //	A += alpha * x * yᴴ + conj(alpha) * y * xᴴ,
   426  //
   427  // where A is an n×n Hermitian matrix, x and y are vectors, and alpha is a scalar.
   428  func Her2(alpha complex64, x, y Vector, a Hermitian) {
   429  	cblas64.Cher2(a.Uplo, a.N, alpha, x.Data, x.Inc, y.Data, y.Inc, a.Data, a.Stride)
   430  }
   431  
   432  // Hpr2 performs a rank-2 update
   433  //
   434  //	A += alpha * x * yᴴ + conj(alpha) * y * xᴴ,
   435  //
   436  // where A is an n×n Hermitian matrix in packed format, x and y are vectors,
   437  // and alpha is a scalar.
   438  func Hpr2(alpha complex64, x, y Vector, a HermitianPacked) {
   439  	cblas64.Chpr2(a.Uplo, a.N, alpha, x.Data, x.Inc, y.Data, y.Inc, a.Data)
   440  }
   441  
   442  // Level 3
   443  
   444  // Gemm computes
   445  //
   446  //	C = alpha * A * B + beta * C,
   447  //
   448  // where A, B, and C are dense matrices, and alpha and beta are scalars.
   449  // tA and tB specify whether A or B are transposed or conjugated.
   450  func Gemm(tA, tB blas.Transpose, alpha complex64, a, b General, beta complex64, c General) {
   451  	var m, n, k int
   452  	if tA == blas.NoTrans {
   453  		m, k = a.Rows, a.Cols
   454  	} else {
   455  		m, k = a.Cols, a.Rows
   456  	}
   457  	if tB == blas.NoTrans {
   458  		n = b.Cols
   459  	} else {
   460  		n = b.Rows
   461  	}
   462  	cblas64.Cgemm(tA, tB, m, n, k, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   463  }
   464  
   465  // Symm performs
   466  //
   467  //	C = alpha * A * B + beta * C  if s == blas.Left,
   468  //	C = alpha * B * A + beta * C  if s == blas.Right,
   469  //
   470  // where A is an n×n or m×m symmetric matrix, B and C are m×n matrices, and
   471  // alpha and beta are scalars.
   472  func Symm(s blas.Side, alpha complex64, a Symmetric, b General, beta complex64, c General) {
   473  	var m, n int
   474  	if s == blas.Left {
   475  		m, n = a.N, b.Cols
   476  	} else {
   477  		m, n = b.Rows, a.N
   478  	}
   479  	cblas64.Csymm(s, a.Uplo, m, n, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   480  }
   481  
   482  // Syrk performs a symmetric rank-k update
   483  //
   484  //	C = alpha * A * Aᵀ + beta * C  if t == blas.NoTrans,
   485  //	C = alpha * Aᵀ * A + beta * C  if t == blas.Trans,
   486  //
   487  // where C is an n×n symmetric matrix, A is an n×k matrix if t == blas.NoTrans
   488  // and a k×n matrix otherwise, and alpha and beta are scalars.
   489  func Syrk(t blas.Transpose, alpha complex64, a General, beta complex64, c Symmetric) {
   490  	var n, k int
   491  	if t == blas.NoTrans {
   492  		n, k = a.Rows, a.Cols
   493  	} else {
   494  		n, k = a.Cols, a.Rows
   495  	}
   496  	cblas64.Csyrk(c.Uplo, t, n, k, alpha, a.Data, a.Stride, beta, c.Data, c.Stride)
   497  }
   498  
   499  // Syr2k performs a symmetric rank-2k update
   500  //
   501  //	C = alpha * A * Bᵀ + alpha * B * Aᵀ + beta * C  if t == blas.NoTrans,
   502  //	C = alpha * Aᵀ * B + alpha * Bᵀ * A + beta * C  if t == blas.Trans,
   503  //
   504  // where C is an n×n symmetric matrix, A and B are n×k matrices if
   505  // t == blas.NoTrans and k×n otherwise, and alpha and beta are scalars.
   506  func Syr2k(t blas.Transpose, alpha complex64, a, b General, beta complex64, c Symmetric) {
   507  	var n, k int
   508  	if t == blas.NoTrans {
   509  		n, k = a.Rows, a.Cols
   510  	} else {
   511  		n, k = a.Cols, a.Rows
   512  	}
   513  	cblas64.Csyr2k(c.Uplo, t, n, k, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   514  }
   515  
   516  // Trmm performs
   517  //
   518  //	B = alpha * A * B   if tA == blas.NoTrans and s == blas.Left,
   519  //	B = alpha * Aᵀ * B  if tA == blas.Trans and s == blas.Left,
   520  //	B = alpha * Aᴴ * B  if tA == blas.ConjTrans and s == blas.Left,
   521  //	B = alpha * B * A   if tA == blas.NoTrans and s == blas.Right,
   522  //	B = alpha * B * Aᵀ  if tA == blas.Trans and s == blas.Right,
   523  //	B = alpha * B * Aᴴ  if tA == blas.ConjTrans and s == blas.Right,
   524  //
   525  // where A is an n×n or m×m triangular matrix, B is an m×n matrix, and alpha is
   526  // a scalar.
   527  func Trmm(s blas.Side, tA blas.Transpose, alpha complex64, a Triangular, b General) {
   528  	cblas64.Ctrmm(s, a.Uplo, tA, a.Diag, b.Rows, b.Cols, alpha, a.Data, a.Stride, b.Data, b.Stride)
   529  }
   530  
   531  // Trsm solves
   532  //
   533  //	A * X = alpha * B   if tA == blas.NoTrans and s == blas.Left,
   534  //	Aᵀ * X = alpha * B  if tA == blas.Trans and s == blas.Left,
   535  //	Aᴴ * X = alpha * B  if tA == blas.ConjTrans and s == blas.Left,
   536  //	X * A = alpha * B   if tA == blas.NoTrans and s == blas.Right,
   537  //	X * Aᵀ = alpha * B  if tA == blas.Trans and s == blas.Right,
   538  //	X * Aᴴ = alpha * B  if tA == blas.ConjTrans and s == blas.Right,
   539  //
   540  // where A is an n×n or m×m triangular matrix, X and B are m×n matrices, and
   541  // alpha is a scalar.
   542  //
   543  // At entry to the function, b contains the values of B, and the result is
   544  // stored in-place into b.
   545  //
   546  // No check is made that A is invertible.
   547  func Trsm(s blas.Side, tA blas.Transpose, alpha complex64, a Triangular, b General) {
   548  	cblas64.Ctrsm(s, a.Uplo, tA, a.Diag, b.Rows, b.Cols, alpha, a.Data, a.Stride, b.Data, b.Stride)
   549  }
   550  
   551  // Hemm performs
   552  //
   553  //	C = alpha * A * B + beta * C  if s == blas.Left,
   554  //	C = alpha * B * A + beta * C  if s == blas.Right,
   555  //
   556  // where A is an n×n or m×m Hermitian matrix, B and C are m×n matrices, and
   557  // alpha and beta are scalars.
   558  func Hemm(s blas.Side, alpha complex64, a Hermitian, b General, beta complex64, c General) {
   559  	var m, n int
   560  	if s == blas.Left {
   561  		m, n = a.N, b.Cols
   562  	} else {
   563  		m, n = b.Rows, a.N
   564  	}
   565  	cblas64.Chemm(s, a.Uplo, m, n, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   566  }
   567  
   568  // Herk performs the Hermitian rank-k update
   569  //
   570  //	C = alpha * A * Aᴴ + beta*C  if t == blas.NoTrans,
   571  //	C = alpha * Aᴴ * A + beta*C  if t == blas.ConjTrans,
   572  //
   573  // where C is an n×n Hermitian matrix, A is an n×k matrix if t == blas.NoTrans
   574  // and a k×n matrix otherwise, and alpha and beta are scalars.
   575  func Herk(t blas.Transpose, alpha float32, a General, beta float32, c Hermitian) {
   576  	var n, k int
   577  	if t == blas.NoTrans {
   578  		n, k = a.Rows, a.Cols
   579  	} else {
   580  		n, k = a.Cols, a.Rows
   581  	}
   582  	cblas64.Cherk(c.Uplo, t, n, k, alpha, a.Data, a.Stride, beta, c.Data, c.Stride)
   583  }
   584  
   585  // Her2k performs the Hermitian rank-2k update
   586  //
   587  //	C = alpha * A * Bᴴ + conj(alpha) * B * Aᴴ + beta * C  if t == blas.NoTrans,
   588  //	C = alpha * Aᴴ * B + conj(alpha) * Bᴴ * A + beta * C  if t == blas.ConjTrans,
   589  //
   590  // where C is an n×n Hermitian matrix, A and B are n×k matrices if t == NoTrans
   591  // and k×n matrices otherwise, and alpha and beta are scalars.
   592  func Her2k(t blas.Transpose, alpha complex64, a, b General, beta float32, c Hermitian) {
   593  	var n, k int
   594  	if t == blas.NoTrans {
   595  		n, k = a.Rows, a.Cols
   596  	} else {
   597  		n, k = a.Cols, a.Rows
   598  	}
   599  	cblas64.Cher2k(c.Uplo, t, n, k, alpha, a.Data, a.Stride, b.Data, b.Stride, beta, c.Data, c.Stride)
   600  }