gonum.org/v1/gonum@v0.15.1-0.20240517103525-f853624cb1bb/mat/matrix.go (about)

     1  // Copyright ©2013 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 mat
     6  
     7  import (
     8  	"math"
     9  
    10  	"gonum.org/v1/gonum/blas"
    11  	"gonum.org/v1/gonum/blas/blas64"
    12  	"gonum.org/v1/gonum/floats/scalar"
    13  	"gonum.org/v1/gonum/lapack"
    14  )
    15  
    16  // Matrix is the basic matrix interface type.
    17  type Matrix interface {
    18  	// Dims returns the dimensions of a Matrix.
    19  	Dims() (r, c int)
    20  
    21  	// At returns the value of a matrix element at row i, column j.
    22  	// It will panic if i or j are out of bounds for the matrix.
    23  	At(i, j int) float64
    24  
    25  	// T returns the transpose of the Matrix. Whether T returns a copy of the
    26  	// underlying data is implementation dependent.
    27  	// This method may be implemented using the Transpose type, which
    28  	// provides an implicit matrix transpose.
    29  	T() Matrix
    30  }
    31  
    32  // allMatrix represents the extra set of methods that all mat Matrix types
    33  // should satisfy. This is used to enforce compile-time consistency between the
    34  // Dense types, especially helpful when adding new features.
    35  type allMatrix interface {
    36  	Reseter
    37  	IsEmpty() bool
    38  	Zero()
    39  }
    40  
    41  // denseMatrix represents the extra set of methods that all Dense Matrix types
    42  // should satisfy. This is used to enforce compile-time consistency between the
    43  // Dense types, especially helpful when adding new features.
    44  type denseMatrix interface {
    45  	DiagView() Diagonal
    46  	Tracer
    47  	Normer
    48  }
    49  
    50  var (
    51  	_ Matrix       = Transpose{}
    52  	_ Untransposer = Transpose{}
    53  )
    54  
    55  // Transpose is a type for performing an implicit matrix transpose. It implements
    56  // the Matrix interface, returning values from the transpose of the matrix within.
    57  type Transpose struct {
    58  	Matrix Matrix
    59  }
    60  
    61  // At returns the value of the element at row i and column j of the transposed
    62  // matrix, that is, row j and column i of the Matrix field.
    63  func (t Transpose) At(i, j int) float64 {
    64  	return t.Matrix.At(j, i)
    65  }
    66  
    67  // Dims returns the dimensions of the transposed matrix. The number of rows returned
    68  // is the number of columns in the Matrix field, and the number of columns is
    69  // the number of rows in the Matrix field.
    70  func (t Transpose) Dims() (r, c int) {
    71  	c, r = t.Matrix.Dims()
    72  	return r, c
    73  }
    74  
    75  // T performs an implicit transpose by returning the Matrix field.
    76  func (t Transpose) T() Matrix {
    77  	return t.Matrix
    78  }
    79  
    80  // Untranspose returns the Matrix field.
    81  func (t Transpose) Untranspose() Matrix {
    82  	return t.Matrix
    83  }
    84  
    85  // Untransposer is a type that can undo an implicit transpose.
    86  type Untransposer interface {
    87  	// Note: This interface is needed to unify all of the Transpose types. In
    88  	// the mat methods, we need to test if the Matrix has been implicitly
    89  	// transposed. If this is checked by testing for the specific Transpose type
    90  	// then the behavior will be different if the user uses T() or TTri() for a
    91  	// triangular matrix.
    92  
    93  	// Untranspose returns the underlying Matrix stored for the implicit transpose.
    94  	Untranspose() Matrix
    95  }
    96  
    97  // UntransposeBander is a type that can undo an implicit band transpose.
    98  type UntransposeBander interface {
    99  	// Untranspose returns the underlying Banded stored for the implicit transpose.
   100  	UntransposeBand() Banded
   101  }
   102  
   103  // UntransposeTrier is a type that can undo an implicit triangular transpose.
   104  type UntransposeTrier interface {
   105  	// Untranspose returns the underlying Triangular stored for the implicit transpose.
   106  	UntransposeTri() Triangular
   107  }
   108  
   109  // UntransposeTriBander is a type that can undo an implicit triangular banded
   110  // transpose.
   111  type UntransposeTriBander interface {
   112  	// Untranspose returns the underlying Triangular stored for the implicit transpose.
   113  	UntransposeTriBand() TriBanded
   114  }
   115  
   116  // Mutable is a matrix interface type that allows elements to be altered.
   117  type Mutable interface {
   118  	// Set alters the matrix element at row i, column j to v.
   119  	// It will panic if i or j are out of bounds for the matrix.
   120  	Set(i, j int, v float64)
   121  
   122  	Matrix
   123  }
   124  
   125  // A RowViewer can return a Vector reflecting a row that is backed by the matrix
   126  // data. The Vector returned will have length equal to the number of columns.
   127  type RowViewer interface {
   128  	RowView(i int) Vector
   129  }
   130  
   131  // A RawRowViewer can return a slice of float64 reflecting a row that is backed by the matrix
   132  // data.
   133  type RawRowViewer interface {
   134  	RawRowView(i int) []float64
   135  }
   136  
   137  // A ColViewer can return a Vector reflecting a column that is backed by the matrix
   138  // data. The Vector returned will have length equal to the number of rows.
   139  type ColViewer interface {
   140  	ColView(j int) Vector
   141  }
   142  
   143  // A RawColViewer can return a slice of float64 reflecting a column that is backed by the matrix
   144  // data.
   145  type RawColViewer interface {
   146  	RawColView(j int) []float64
   147  }
   148  
   149  // A ClonerFrom can make a copy of a into the receiver, overwriting the previous value of the
   150  // receiver. The clone operation does not make any restriction on shape and will not cause
   151  // shadowing.
   152  type ClonerFrom interface {
   153  	CloneFrom(a Matrix)
   154  }
   155  
   156  // A Reseter can reset the matrix so that it can be reused as the receiver of a dimensionally
   157  // restricted operation. This is commonly used when the matrix is being used as a workspace
   158  // or temporary matrix.
   159  //
   160  // If the matrix is a view, using Reset may result in data corruption in elements outside
   161  // the view. Similarly, if the matrix shares backing data with another variable, using
   162  // Reset may lead to unexpected changes in data values.
   163  type Reseter interface {
   164  	Reset()
   165  }
   166  
   167  // A Copier can make a copy of elements of a into the receiver. The submatrix copied
   168  // starts at row and column 0 and has dimensions equal to the minimum dimensions of
   169  // the two matrices. The number of row and columns copied is returned.
   170  // Copy will copy from a source that aliases the receiver unless the source is transposed;
   171  // an aliasing transpose copy will panic with the exception for a special case when
   172  // the source data has a unitary increment or stride.
   173  type Copier interface {
   174  	Copy(a Matrix) (r, c int)
   175  }
   176  
   177  // A Grower can grow the size of the represented matrix by the given number of rows and columns.
   178  // Growing beyond the size given by the Caps method will result in the allocation of a new
   179  // matrix and copying of the elements. If Grow is called with negative increments it will
   180  // panic with ErrIndexOutOfRange.
   181  type Grower interface {
   182  	Caps() (r, c int)
   183  	Grow(r, c int) Matrix
   184  }
   185  
   186  // A RawMatrixSetter can set the underlying blas64.General used by the receiver. There is no restriction
   187  // on the shape of the receiver. Changes to the receiver's elements will be reflected in the blas64.General.Data.
   188  type RawMatrixSetter interface {
   189  	SetRawMatrix(a blas64.General)
   190  }
   191  
   192  // A RawMatrixer can return a blas64.General representation of the receiver. Changes to the blas64.General.Data
   193  // slice will be reflected in the original matrix, changes to the Rows, Cols and Stride fields will not.
   194  type RawMatrixer interface {
   195  	RawMatrix() blas64.General
   196  }
   197  
   198  // A RawVectorer can return a blas64.Vector representation of the receiver. Changes to the blas64.Vector.Data
   199  // slice will be reflected in the original matrix, changes to the Inc field will not.
   200  type RawVectorer interface {
   201  	RawVector() blas64.Vector
   202  }
   203  
   204  // A NonZeroDoer can call a function for each non-zero element of the receiver.
   205  // The parameters of the function are the element indices and its value.
   206  type NonZeroDoer interface {
   207  	DoNonZero(func(i, j int, v float64))
   208  }
   209  
   210  // A RowNonZeroDoer can call a function for each non-zero element of a row of the receiver.
   211  // The parameters of the function are the element indices and its value.
   212  type RowNonZeroDoer interface {
   213  	DoRowNonZero(i int, fn func(i, j int, v float64))
   214  }
   215  
   216  // A ColNonZeroDoer can call a function for each non-zero element of a column of the receiver.
   217  // The parameters of the function are the element indices and its value.
   218  type ColNonZeroDoer interface {
   219  	DoColNonZero(j int, fn func(i, j int, v float64))
   220  }
   221  
   222  // A SolveToer can solve a linear system A⋅X = B or Aᵀ⋅X = B where A is a matrix
   223  // represented by the receiver and B is a given matrix, storing the result into
   224  // dst.
   225  //
   226  // If dst is empty, SolveTo will resize it to the correct size, otherwise it
   227  // must have the correct size. Individual implementations may impose other
   228  // restrictions on the input parameters, for example that A is a square matrix.
   229  type SolveToer interface {
   230  	SolveTo(dst *Dense, trans bool, b Matrix) error
   231  }
   232  
   233  // untranspose untransposes a matrix if applicable. If a is an Untransposer, then
   234  // untranspose returns the underlying matrix and true. If it is not, then it returns
   235  // the input matrix and false.
   236  func untranspose(a Matrix) (Matrix, bool) {
   237  	if ut, ok := a.(Untransposer); ok {
   238  		return ut.Untranspose(), true
   239  	}
   240  	return a, false
   241  }
   242  
   243  // untransposeExtract returns an untransposed matrix in a built-in matrix type.
   244  //
   245  // The untransposed matrix is returned unaltered if it is a built-in matrix type.
   246  // Otherwise, if it implements a Raw method, an appropriate built-in type value
   247  // is returned holding the raw matrix value of the input. If neither of these
   248  // is possible, the untransposed matrix is returned.
   249  func untransposeExtract(a Matrix) (Matrix, bool) {
   250  	ut, trans := untranspose(a)
   251  	switch m := ut.(type) {
   252  	case *DiagDense, *SymBandDense, *TriBandDense, *BandDense, *TriDense, *SymDense, *Dense, *VecDense, *Tridiag:
   253  		return m, trans
   254  	// TODO(btracey): Add here if we ever have an equivalent of RawDiagDense.
   255  	case RawSymBander:
   256  		rsb := m.RawSymBand()
   257  		if rsb.Uplo != blas.Upper {
   258  			return ut, trans
   259  		}
   260  		var sb SymBandDense
   261  		sb.SetRawSymBand(rsb)
   262  		return &sb, trans
   263  	case RawTriBander:
   264  		rtb := m.RawTriBand()
   265  		if rtb.Diag == blas.Unit {
   266  			return ut, trans
   267  		}
   268  		var tb TriBandDense
   269  		tb.SetRawTriBand(rtb)
   270  		return &tb, trans
   271  	case RawBander:
   272  		var b BandDense
   273  		b.SetRawBand(m.RawBand())
   274  		return &b, trans
   275  	case RawTriangular:
   276  		rt := m.RawTriangular()
   277  		if rt.Diag == blas.Unit {
   278  			return ut, trans
   279  		}
   280  		var t TriDense
   281  		t.SetRawTriangular(rt)
   282  		return &t, trans
   283  	case RawSymmetricer:
   284  		rs := m.RawSymmetric()
   285  		if rs.Uplo != blas.Upper {
   286  			return ut, trans
   287  		}
   288  		var s SymDense
   289  		s.SetRawSymmetric(rs)
   290  		return &s, trans
   291  	case RawMatrixer:
   292  		var d Dense
   293  		d.SetRawMatrix(m.RawMatrix())
   294  		return &d, trans
   295  	case RawVectorer:
   296  		var v VecDense
   297  		v.SetRawVector(m.RawVector())
   298  		return &v, trans
   299  	case RawTridiagonaler:
   300  		var d Tridiag
   301  		d.SetRawTridiagonal(m.RawTridiagonal())
   302  		return &d, trans
   303  	default:
   304  		return ut, trans
   305  	}
   306  }
   307  
   308  // TODO(btracey): Consider adding CopyCol/CopyRow if the behavior seems useful.
   309  // TODO(btracey): Add in fast paths to Row/Col for the other concrete types
   310  // (TriDense, etc.) as well as relevant interfaces (RowColer, RawRowViewer, etc.)
   311  
   312  // Col copies the elements in the jth column of the matrix into the slice dst.
   313  // The length of the provided slice must equal the number of rows, unless the
   314  // slice is nil in which case a new slice is first allocated.
   315  func Col(dst []float64, j int, a Matrix) []float64 {
   316  	r, c := a.Dims()
   317  	if j < 0 || j >= c {
   318  		panic(ErrColAccess)
   319  	}
   320  	if dst == nil {
   321  		dst = make([]float64, r)
   322  	} else {
   323  		if len(dst) != r {
   324  			panic(ErrColLength)
   325  		}
   326  	}
   327  	aU, aTrans := untranspose(a)
   328  	if rm, ok := aU.(RawMatrixer); ok {
   329  		m := rm.RawMatrix()
   330  		if aTrans {
   331  			copy(dst, m.Data[j*m.Stride:j*m.Stride+m.Cols])
   332  			return dst
   333  		}
   334  		blas64.Copy(blas64.Vector{N: r, Inc: m.Stride, Data: m.Data[j:]},
   335  			blas64.Vector{N: r, Inc: 1, Data: dst},
   336  		)
   337  		return dst
   338  	}
   339  	for i := 0; i < r; i++ {
   340  		dst[i] = a.At(i, j)
   341  	}
   342  	return dst
   343  }
   344  
   345  // Row copies the elements in the ith row of the matrix into the slice dst.
   346  // The length of the provided slice must equal the number of columns, unless the
   347  // slice is nil in which case a new slice is first allocated.
   348  func Row(dst []float64, i int, a Matrix) []float64 {
   349  	r, c := a.Dims()
   350  	if i < 0 || i >= r {
   351  		panic(ErrColAccess)
   352  	}
   353  	if dst == nil {
   354  		dst = make([]float64, c)
   355  	} else {
   356  		if len(dst) != c {
   357  			panic(ErrRowLength)
   358  		}
   359  	}
   360  	aU, aTrans := untranspose(a)
   361  	if rm, ok := aU.(RawMatrixer); ok {
   362  		m := rm.RawMatrix()
   363  		if aTrans {
   364  			blas64.Copy(blas64.Vector{N: c, Inc: m.Stride, Data: m.Data[i:]},
   365  				blas64.Vector{N: c, Inc: 1, Data: dst},
   366  			)
   367  			return dst
   368  		}
   369  		copy(dst, m.Data[i*m.Stride:i*m.Stride+m.Cols])
   370  		return dst
   371  	}
   372  	for j := 0; j < c; j++ {
   373  		dst[j] = a.At(i, j)
   374  	}
   375  	return dst
   376  }
   377  
   378  // Cond returns the condition number of the given matrix under the given norm.
   379  // The condition number must be based on the 1-norm, 2-norm or ∞-norm.
   380  // Cond will panic with ErrZeroLength if the matrix has zero size.
   381  //
   382  // BUG(btracey): The computation of the 1-norm and ∞-norm for non-square matrices
   383  // is inaccurate, although is typically the right order of magnitude. See
   384  // https://github.com/xianyi/OpenBLAS/issues/636. While the value returned will
   385  // change with the resolution of this bug, the result from Cond will match the
   386  // condition number used internally.
   387  func Cond(a Matrix, norm float64) float64 {
   388  	m, n := a.Dims()
   389  	if m == 0 || n == 0 {
   390  		panic(ErrZeroLength)
   391  	}
   392  	var lnorm lapack.MatrixNorm
   393  	switch norm {
   394  	default:
   395  		panic("mat: bad norm value")
   396  	case 1:
   397  		lnorm = lapack.MaxColumnSum
   398  	case 2:
   399  		var svd SVD
   400  		ok := svd.Factorize(a, SVDNone)
   401  		if !ok {
   402  			return math.Inf(1)
   403  		}
   404  		return svd.Cond()
   405  	case math.Inf(1):
   406  		lnorm = lapack.MaxRowSum
   407  	}
   408  
   409  	if m == n {
   410  		// Use the LU decomposition to compute the condition number.
   411  		var lu LU
   412  		lu.factorize(a, lnorm)
   413  		return lu.Cond()
   414  	}
   415  	if m > n {
   416  		// Use the QR factorization to compute the condition number.
   417  		var qr QR
   418  		qr.factorize(a, lnorm)
   419  		return qr.Cond()
   420  	}
   421  	// Use the LQ factorization to compute the condition number.
   422  	var lq LQ
   423  	lq.factorize(a, lnorm)
   424  	return lq.Cond()
   425  }
   426  
   427  // Det returns the determinant of the square matrix a. In many expressions using
   428  // LogDet will be more numerically stable.
   429  //
   430  // Det panics with ErrSquare if a is not square and with ErrZeroLength if a has
   431  // zero size.
   432  func Det(a Matrix) float64 {
   433  	det, sign := LogDet(a)
   434  	return math.Exp(det) * sign
   435  }
   436  
   437  // Dot returns the sum of the element-wise product of a and b.
   438  //
   439  // Dot panics with ErrShape if the vector sizes are unequal and with
   440  // ErrZeroLength if the sizes are zero.
   441  func Dot(a, b Vector) float64 {
   442  	la := a.Len()
   443  	lb := b.Len()
   444  	if la != lb {
   445  		panic(ErrShape)
   446  	}
   447  	if la == 0 {
   448  		panic(ErrZeroLength)
   449  	}
   450  	if arv, ok := a.(RawVectorer); ok {
   451  		if brv, ok := b.(RawVectorer); ok {
   452  			return blas64.Dot(arv.RawVector(), brv.RawVector())
   453  		}
   454  	}
   455  	var sum float64
   456  	for i := 0; i < la; i++ {
   457  		sum += a.At(i, 0) * b.At(i, 0)
   458  	}
   459  	return sum
   460  }
   461  
   462  // Equal returns whether the matrices a and b have the same size
   463  // and are element-wise equal.
   464  func Equal(a, b Matrix) bool {
   465  	ar, ac := a.Dims()
   466  	br, bc := b.Dims()
   467  	if ar != br || ac != bc {
   468  		return false
   469  	}
   470  	aU, aTrans := untranspose(a)
   471  	bU, bTrans := untranspose(b)
   472  	if rma, ok := aU.(RawMatrixer); ok {
   473  		if rmb, ok := bU.(RawMatrixer); ok {
   474  			ra := rma.RawMatrix()
   475  			rb := rmb.RawMatrix()
   476  			if aTrans == bTrans {
   477  				for i := 0; i < ra.Rows; i++ {
   478  					for j := 0; j < ra.Cols; j++ {
   479  						if ra.Data[i*ra.Stride+j] != rb.Data[i*rb.Stride+j] {
   480  							return false
   481  						}
   482  					}
   483  				}
   484  				return true
   485  			}
   486  			for i := 0; i < ra.Rows; i++ {
   487  				for j := 0; j < ra.Cols; j++ {
   488  					if ra.Data[i*ra.Stride+j] != rb.Data[j*rb.Stride+i] {
   489  						return false
   490  					}
   491  				}
   492  			}
   493  			return true
   494  		}
   495  	}
   496  	if rma, ok := aU.(RawSymmetricer); ok {
   497  		if rmb, ok := bU.(RawSymmetricer); ok {
   498  			ra := rma.RawSymmetric()
   499  			rb := rmb.RawSymmetric()
   500  			// Symmetric matrices are always upper and equal to their transpose.
   501  			for i := 0; i < ra.N; i++ {
   502  				for j := i; j < ra.N; j++ {
   503  					if ra.Data[i*ra.Stride+j] != rb.Data[i*rb.Stride+j] {
   504  						return false
   505  					}
   506  				}
   507  			}
   508  			return true
   509  		}
   510  	}
   511  	if ra, ok := aU.(*VecDense); ok {
   512  		if rb, ok := bU.(*VecDense); ok {
   513  			// If the raw vectors are the same length they must either both be
   514  			// transposed or both not transposed (or have length 1).
   515  			for i := 0; i < ra.mat.N; i++ {
   516  				if ra.mat.Data[i*ra.mat.Inc] != rb.mat.Data[i*rb.mat.Inc] {
   517  					return false
   518  				}
   519  			}
   520  			return true
   521  		}
   522  	}
   523  	for i := 0; i < ar; i++ {
   524  		for j := 0; j < ac; j++ {
   525  			if a.At(i, j) != b.At(i, j) {
   526  				return false
   527  			}
   528  		}
   529  	}
   530  	return true
   531  }
   532  
   533  // EqualApprox returns whether the matrices a and b have the same size and contain all equal
   534  // elements with tolerance for element-wise equality specified by epsilon. Matrices
   535  // with non-equal shapes are not equal.
   536  func EqualApprox(a, b Matrix, epsilon float64) bool {
   537  	ar, ac := a.Dims()
   538  	br, bc := b.Dims()
   539  	if ar != br || ac != bc {
   540  		return false
   541  	}
   542  	aU, aTrans := untranspose(a)
   543  	bU, bTrans := untranspose(b)
   544  	if rma, ok := aU.(RawMatrixer); ok {
   545  		if rmb, ok := bU.(RawMatrixer); ok {
   546  			ra := rma.RawMatrix()
   547  			rb := rmb.RawMatrix()
   548  			if aTrans == bTrans {
   549  				for i := 0; i < ra.Rows; i++ {
   550  					for j := 0; j < ra.Cols; j++ {
   551  						if !scalar.EqualWithinAbsOrRel(ra.Data[i*ra.Stride+j], rb.Data[i*rb.Stride+j], epsilon, epsilon) {
   552  							return false
   553  						}
   554  					}
   555  				}
   556  				return true
   557  			}
   558  			for i := 0; i < ra.Rows; i++ {
   559  				for j := 0; j < ra.Cols; j++ {
   560  					if !scalar.EqualWithinAbsOrRel(ra.Data[i*ra.Stride+j], rb.Data[j*rb.Stride+i], epsilon, epsilon) {
   561  						return false
   562  					}
   563  				}
   564  			}
   565  			return true
   566  		}
   567  	}
   568  	if rma, ok := aU.(RawSymmetricer); ok {
   569  		if rmb, ok := bU.(RawSymmetricer); ok {
   570  			ra := rma.RawSymmetric()
   571  			rb := rmb.RawSymmetric()
   572  			// Symmetric matrices are always upper and equal to their transpose.
   573  			for i := 0; i < ra.N; i++ {
   574  				for j := i; j < ra.N; j++ {
   575  					if !scalar.EqualWithinAbsOrRel(ra.Data[i*ra.Stride+j], rb.Data[i*rb.Stride+j], epsilon, epsilon) {
   576  						return false
   577  					}
   578  				}
   579  			}
   580  			return true
   581  		}
   582  	}
   583  	if ra, ok := aU.(*VecDense); ok {
   584  		if rb, ok := bU.(*VecDense); ok {
   585  			// If the raw vectors are the same length they must either both be
   586  			// transposed or both not transposed (or have length 1).
   587  			for i := 0; i < ra.mat.N; i++ {
   588  				if !scalar.EqualWithinAbsOrRel(ra.mat.Data[i*ra.mat.Inc], rb.mat.Data[i*rb.mat.Inc], epsilon, epsilon) {
   589  					return false
   590  				}
   591  			}
   592  			return true
   593  		}
   594  	}
   595  	for i := 0; i < ar; i++ {
   596  		for j := 0; j < ac; j++ {
   597  			if !scalar.EqualWithinAbsOrRel(a.At(i, j), b.At(i, j), epsilon, epsilon) {
   598  				return false
   599  			}
   600  		}
   601  	}
   602  	return true
   603  }
   604  
   605  // LogDet returns the log of the determinant and the sign of the determinant
   606  // for the matrix that has been factorized. Numerical stability in product and
   607  // division expressions is generally improved by working in log space.
   608  //
   609  // LogDet panics with ErrSquare is a is not square and with ErrZeroLength if a
   610  // has zero size.
   611  func LogDet(a Matrix) (det float64, sign float64) {
   612  	// TODO(btracey): Add specialized routines for TriDense, etc.
   613  	var lu LU
   614  	lu.Factorize(a)
   615  	return lu.LogDet()
   616  }
   617  
   618  // Max returns the largest element value of the matrix A.
   619  //
   620  // Max will panic with ErrZeroLength if the matrix has zero size.
   621  func Max(a Matrix) float64 {
   622  	r, c := a.Dims()
   623  	if r == 0 || c == 0 {
   624  		panic(ErrZeroLength)
   625  	}
   626  	// Max(A) = Max(Aᵀ)
   627  	aU, _ := untranspose(a)
   628  	switch m := aU.(type) {
   629  	case RawMatrixer:
   630  		rm := m.RawMatrix()
   631  		max := math.Inf(-1)
   632  		for i := 0; i < rm.Rows; i++ {
   633  			for _, v := range rm.Data[i*rm.Stride : i*rm.Stride+rm.Cols] {
   634  				if v > max {
   635  					max = v
   636  				}
   637  			}
   638  		}
   639  		return max
   640  	case RawTriangular:
   641  		rm := m.RawTriangular()
   642  		// The max of a triangular is at least 0 unless the size is 1.
   643  		if rm.N == 1 {
   644  			return rm.Data[0]
   645  		}
   646  		max := 0.0
   647  		if rm.Uplo == blas.Upper {
   648  			for i := 0; i < rm.N; i++ {
   649  				for _, v := range rm.Data[i*rm.Stride+i : i*rm.Stride+rm.N] {
   650  					if v > max {
   651  						max = v
   652  					}
   653  				}
   654  			}
   655  			return max
   656  		}
   657  		for i := 0; i < rm.N; i++ {
   658  			for _, v := range rm.Data[i*rm.Stride : i*rm.Stride+i+1] {
   659  				if v > max {
   660  					max = v
   661  				}
   662  			}
   663  		}
   664  		return max
   665  	case RawSymmetricer:
   666  		rm := m.RawSymmetric()
   667  		if rm.Uplo != blas.Upper {
   668  			panic(badSymTriangle)
   669  		}
   670  		max := math.Inf(-1)
   671  		for i := 0; i < rm.N; i++ {
   672  			for _, v := range rm.Data[i*rm.Stride+i : i*rm.Stride+rm.N] {
   673  				if v > max {
   674  					max = v
   675  				}
   676  			}
   677  		}
   678  		return max
   679  	default:
   680  		r, c := aU.Dims()
   681  		max := math.Inf(-1)
   682  		for i := 0; i < r; i++ {
   683  			for j := 0; j < c; j++ {
   684  				v := aU.At(i, j)
   685  				if v > max {
   686  					max = v
   687  				}
   688  			}
   689  		}
   690  		return max
   691  	}
   692  }
   693  
   694  // Min returns the smallest element value of the matrix A.
   695  //
   696  // Min will panic with ErrZeroLength if the matrix has zero size.
   697  func Min(a Matrix) float64 {
   698  	r, c := a.Dims()
   699  	if r == 0 || c == 0 {
   700  		panic(ErrZeroLength)
   701  	}
   702  	// Min(A) = Min(Aᵀ)
   703  	aU, _ := untranspose(a)
   704  	switch m := aU.(type) {
   705  	case RawMatrixer:
   706  		rm := m.RawMatrix()
   707  		min := math.Inf(1)
   708  		for i := 0; i < rm.Rows; i++ {
   709  			for _, v := range rm.Data[i*rm.Stride : i*rm.Stride+rm.Cols] {
   710  				if v < min {
   711  					min = v
   712  				}
   713  			}
   714  		}
   715  		return min
   716  	case RawTriangular:
   717  		rm := m.RawTriangular()
   718  		// The min of a triangular is at most 0 unless the size is 1.
   719  		if rm.N == 1 {
   720  			return rm.Data[0]
   721  		}
   722  		min := 0.0
   723  		if rm.Uplo == blas.Upper {
   724  			for i := 0; i < rm.N; i++ {
   725  				for _, v := range rm.Data[i*rm.Stride+i : i*rm.Stride+rm.N] {
   726  					if v < min {
   727  						min = v
   728  					}
   729  				}
   730  			}
   731  			return min
   732  		}
   733  		for i := 0; i < rm.N; i++ {
   734  			for _, v := range rm.Data[i*rm.Stride : i*rm.Stride+i+1] {
   735  				if v < min {
   736  					min = v
   737  				}
   738  			}
   739  		}
   740  		return min
   741  	case RawSymmetricer:
   742  		rm := m.RawSymmetric()
   743  		if rm.Uplo != blas.Upper {
   744  			panic(badSymTriangle)
   745  		}
   746  		min := math.Inf(1)
   747  		for i := 0; i < rm.N; i++ {
   748  			for _, v := range rm.Data[i*rm.Stride+i : i*rm.Stride+rm.N] {
   749  				if v < min {
   750  					min = v
   751  				}
   752  			}
   753  		}
   754  		return min
   755  	default:
   756  		r, c := aU.Dims()
   757  		min := math.Inf(1)
   758  		for i := 0; i < r; i++ {
   759  			for j := 0; j < c; j++ {
   760  				v := aU.At(i, j)
   761  				if v < min {
   762  					min = v
   763  				}
   764  			}
   765  		}
   766  		return min
   767  	}
   768  }
   769  
   770  // A Normer can compute a norm of the matrix. Valid norms are:
   771  //
   772  //	1 - The maximum absolute column sum
   773  //	2 - The Frobenius norm, the square root of the sum of the squares of the elements
   774  //	Inf - The maximum absolute row sum
   775  type Normer interface {
   776  	Norm(norm float64) float64
   777  }
   778  
   779  // Norm returns the specified norm of the matrix A. Valid norms are:
   780  //
   781  //	1 - The maximum absolute column sum
   782  //	2 - The Frobenius norm, the square root of the sum of the squares of the elements
   783  //	Inf - The maximum absolute row sum
   784  //
   785  // If a is a Normer, its Norm method will be used to calculate the norm.
   786  //
   787  // Norm will panic with ErrNormOrder if an illegal norm is specified and with
   788  // ErrShape if the matrix has zero size.
   789  func Norm(a Matrix, norm float64) float64 {
   790  	r, c := a.Dims()
   791  	if r == 0 || c == 0 {
   792  		panic(ErrZeroLength)
   793  	}
   794  	m, trans := untransposeExtract(a)
   795  	if m, ok := m.(Normer); ok {
   796  		if trans {
   797  			switch norm {
   798  			case 1:
   799  				norm = math.Inf(1)
   800  			case math.Inf(1):
   801  				norm = 1
   802  			}
   803  		}
   804  		return m.Norm(norm)
   805  	}
   806  	switch norm {
   807  	default:
   808  		panic(ErrNormOrder)
   809  	case 1:
   810  		var max float64
   811  		for j := 0; j < c; j++ {
   812  			var sum float64
   813  			for i := 0; i < r; i++ {
   814  				sum += math.Abs(a.At(i, j))
   815  			}
   816  			if sum > max {
   817  				max = sum
   818  			}
   819  		}
   820  		return max
   821  	case 2:
   822  		var sum float64
   823  		for i := 0; i < r; i++ {
   824  			for j := 0; j < c; j++ {
   825  				v := a.At(i, j)
   826  				sum += v * v
   827  			}
   828  		}
   829  		return math.Sqrt(sum)
   830  	case math.Inf(1):
   831  		var max float64
   832  		for i := 0; i < r; i++ {
   833  			var sum float64
   834  			for j := 0; j < c; j++ {
   835  				sum += math.Abs(a.At(i, j))
   836  			}
   837  			if sum > max {
   838  				max = sum
   839  			}
   840  		}
   841  		return max
   842  	}
   843  }
   844  
   845  // normLapack converts the float64 norm input in Norm to a lapack.MatrixNorm.
   846  func normLapack(norm float64, aTrans bool) lapack.MatrixNorm {
   847  	switch norm {
   848  	case 1:
   849  		n := lapack.MaxColumnSum
   850  		if aTrans {
   851  			n = lapack.MaxRowSum
   852  		}
   853  		return n
   854  	case 2:
   855  		return lapack.Frobenius
   856  	case math.Inf(1):
   857  		n := lapack.MaxRowSum
   858  		if aTrans {
   859  			n = lapack.MaxColumnSum
   860  		}
   861  		return n
   862  	default:
   863  		panic(ErrNormOrder)
   864  	}
   865  }
   866  
   867  // Sum returns the sum of the elements of the matrix.
   868  //
   869  // Sum will panic with ErrZeroLength if the matrix has zero size.
   870  func Sum(a Matrix) float64 {
   871  	r, c := a.Dims()
   872  	if r == 0 || c == 0 {
   873  		panic(ErrZeroLength)
   874  	}
   875  	var sum float64
   876  	aU, _ := untranspose(a)
   877  	switch rma := aU.(type) {
   878  	case RawSymmetricer:
   879  		rm := rma.RawSymmetric()
   880  		for i := 0; i < rm.N; i++ {
   881  			// Diagonals count once while off-diagonals count twice.
   882  			sum += rm.Data[i*rm.Stride+i]
   883  			var s float64
   884  			for _, v := range rm.Data[i*rm.Stride+i+1 : i*rm.Stride+rm.N] {
   885  				s += v
   886  			}
   887  			sum += 2 * s
   888  		}
   889  		return sum
   890  	case RawTriangular:
   891  		rm := rma.RawTriangular()
   892  		var startIdx, endIdx int
   893  		for i := 0; i < rm.N; i++ {
   894  			// Start and end index for this triangle-row.
   895  			switch rm.Uplo {
   896  			case blas.Upper:
   897  				startIdx = i
   898  				endIdx = rm.N
   899  			case blas.Lower:
   900  				startIdx = 0
   901  				endIdx = i + 1
   902  			default:
   903  				panic(badTriangle)
   904  			}
   905  			for _, v := range rm.Data[i*rm.Stride+startIdx : i*rm.Stride+endIdx] {
   906  				sum += v
   907  			}
   908  		}
   909  		return sum
   910  	case RawMatrixer:
   911  		rm := rma.RawMatrix()
   912  		for i := 0; i < rm.Rows; i++ {
   913  			for _, v := range rm.Data[i*rm.Stride : i*rm.Stride+rm.Cols] {
   914  				sum += v
   915  			}
   916  		}
   917  		return sum
   918  	case *VecDense:
   919  		rm := rma.RawVector()
   920  		for i := 0; i < rm.N; i++ {
   921  			sum += rm.Data[i*rm.Inc]
   922  		}
   923  		return sum
   924  	default:
   925  		r, c := a.Dims()
   926  		for i := 0; i < r; i++ {
   927  			for j := 0; j < c; j++ {
   928  				sum += a.At(i, j)
   929  			}
   930  		}
   931  		return sum
   932  	}
   933  }
   934  
   935  // A Tracer can compute the trace of the matrix. Trace must panic with ErrSquare
   936  // if the matrix is not square.
   937  type Tracer interface {
   938  	Trace() float64
   939  }
   940  
   941  // Trace returns the trace of the matrix. If a is a Tracer, its Trace method
   942  // will be used to calculate the matrix trace.
   943  //
   944  // Trace will panic with ErrSquare if the matrix is not square and with
   945  // ErrZeroLength if the matrix has zero size.
   946  func Trace(a Matrix) float64 {
   947  	r, c := a.Dims()
   948  	if r == 0 || c == 0 {
   949  		panic(ErrZeroLength)
   950  	}
   951  	m, _ := untransposeExtract(a)
   952  	if t, ok := m.(Tracer); ok {
   953  		return t.Trace()
   954  	}
   955  	if r != c {
   956  		panic(ErrSquare)
   957  	}
   958  	var v float64
   959  	for i := 0; i < r; i++ {
   960  		v += a.At(i, i)
   961  	}
   962  	return v
   963  }
   964  
   965  // use returns a float64 slice with l elements, using f if it
   966  // has the necessary capacity, otherwise creating a new slice.
   967  func use(f []float64, l int) []float64 {
   968  	if l <= cap(f) {
   969  		return f[:l]
   970  	}
   971  	return make([]float64, l)
   972  }
   973  
   974  // useZeroed returns a float64 slice with l elements, using f if it
   975  // has the necessary capacity, otherwise creating a new slice. The
   976  // elements of the returned slice are guaranteed to be zero.
   977  func useZeroed(f []float64, l int) []float64 {
   978  	if l <= cap(f) {
   979  		f = f[:l]
   980  		zero(f)
   981  		return f
   982  	}
   983  	return make([]float64, l)
   984  }
   985  
   986  // zero zeros the given slice's elements.
   987  func zero(f []float64) {
   988  	for i := range f {
   989  		f[i] = 0
   990  	}
   991  }
   992  
   993  // useInt returns an int slice with l elements, using i if it
   994  // has the necessary capacity, otherwise creating a new slice.
   995  func useInt(i []int, l int) []int {
   996  	if l <= cap(i) {
   997  		return i[:l]
   998  	}
   999  	return make([]int, l)
  1000  }