-
-
Notifications
You must be signed in to change notification settings - Fork 5.6k
/
Copy pathLinearAlgebra.jl
355 lines (311 loc) · 7.99 KB
/
LinearAlgebra.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
# This file is a part of Julia. License is MIT: https://julialang.org/license
__precompile__(true)
"""
Linear algebra module. Provides array arithmetic,
matrix factorizations and other linear algebra related
functionality.
"""
module LinearAlgebra
import Base: \, /, *, ^, +, -, ==
import Base: USE_BLAS64, abs, acos, acosh, acot, acoth, acsc, acsch, adjoint, asec, asech,
asin, asinh, atan, atanh, axes, big, broadcast, ceil, conj, convert, copy, copyto!, cos,
cosh, cot, coth, csc, csch, eltype, exp, findmax, findmin, fill!, floor, getindex, hcat,
getproperty, imag, inv, isapprox, isone, IndexStyle, kron, length, log, map, ndims,
oneunit, parent, power_by_squaring, print_matrix, promote_rule, real, round, sec, sech,
setindex!, show, similar, sin, sincos, sinh, size, size_to_strides, sqrt, StridedReinterpretArray,
StridedReshapedArray, strides, stride, tan, tanh, transpose, trunc, typed_hcat, vec
using Base: hvcat_fill, iszero, IndexLinear, _length, promote_op, promote_typeof,
@propagate_inbounds, @pure, reduce, typed_vcat
# We use `_length` because of non-1 indices; releases after julia 0.5
# can go back to `length`. `_length(A)` is equivalent to `length(linearindices(A))`.
export
# Modules
LAPACK,
BLAS,
# Types
Adjoint,
Transpose,
SymTridiagonal,
Tridiagonal,
Bidiagonal,
Factorization,
BunchKaufman,
Cholesky,
CholeskyPivoted,
Eigen,
GeneralizedEigen,
GeneralizedSVD,
GeneralizedSchur,
Hessenberg,
LU,
LDLt,
QR,
QRPivoted,
LQ,
Schur,
SVD,
Hermitian,
Symmetric,
LowerTriangular,
UpperTriangular,
Diagonal,
UniformScaling,
# Functions
axpy!,
axpby!,
bkfact,
bkfact!,
chol,
cholfact,
cholfact!,
cond,
condskeel,
copyto!,
copy_transpose!,
cross,
adjoint,
adjoint!,
det,
diag,
diagind,
diagm,
diff,
dot,
eig,
eigfact,
eigfact!,
eigmax,
eigmin,
eigvals,
eigvals!,
eigvecs,
factorize,
givens,
hessfact,
hessfact!,
isdiag,
ishermitian,
isposdef,
isposdef!,
issuccess,
issymmetric,
istril,
istriu,
kron,
ldiv!,
ldltfact!,
ldltfact,
linreg,
logabsdet,
logdet,
lu,
lufact,
lufact!,
lyap,
mul!,
mul1!,
mul2!,
norm,
normalize,
normalize!,
nullspace,
ordschur!,
ordschur,
peakflops,
pinv,
qr,
qrfact!,
qrfact,
lq,
lqfact!,
lqfact,
rank,
rdiv!,
schur,
schurfact!,
schurfact,
svd,
svdfact!,
svdfact,
svdvals!,
svdvals,
sylvester,
trace,
transpose,
transpose!,
transpose_type,
tril,
triu,
tril!,
triu!,
vecdot,
vecnorm,
# Operators
\,
/,
# Constants
I
const BlasFloat = Union{Float64,Float32,ComplexF64,ComplexF32}
const BlasReal = Union{Float64,Float32}
const BlasComplex = Union{ComplexF64,ComplexF32}
if USE_BLAS64
const BlasInt = Int64
else
const BlasInt = Int32
end
# Check that stride of matrix/vector is 1
# Writing like this to avoid splatting penalty when called with multiple arguments,
# see PR 16416
"""
stride1(A) -> Int
Return the distance between successive array elements
in dimension 1 in units of element size.
# Examples
```jldoctest
julia> A = [1,2,3,4]
4-element Array{Int64,1}:
1
2
3
4
julia> LinearAlgebra.stride1(A)
1
julia> B = view(A, 2:2:4)
2-element view(::Array{Int64,1}, 2:2:4) with eltype Int64:
2
4
julia> LinearAlgebra.stride1(B)
2
```
"""
stride1(x) = stride(x,1)
stride1(x::Array) = 1
stride1(x::DenseArray) = stride(x, 1)::Int
@inline chkstride1(A...) = _chkstride1(true, A...)
@noinline _chkstride1(ok::Bool) = ok || error("matrix does not have contiguous columns")
@inline _chkstride1(ok::Bool, A, B...) = _chkstride1(ok & (stride1(A) == 1), B...)
"""
LinearAlgebra.checksquare(A)
Check that a matrix is square, then return its common dimension.
For multiple arguments, return a vector.
# Examples
```jldoctest
julia> A = fill(1, (4,4)); B = fill(1, (5,5));
julia> LinearAlgebra.checksquare(A, B)
2-element Array{Int64,1}:
4
5
```
"""
function checksquare(A)
m,n = size(A)
m == n || throw(DimensionMismatch("matrix is not square: dimensions are $(size(A))"))
m
end
function checksquare(A...)
sizes = Int[]
for a in A
size(a,1)==size(a,2) || throw(DimensionMismatch("matrix is not square: dimensions are $(size(a))"))
push!(sizes, size(a,1))
end
return sizes
end
function char_uplo(uplo::Symbol)
if uplo == :U
'U'
elseif uplo == :L
'L'
else
throw(ArgumentError("uplo argument must be either :U (upper) or :L (lower)"))
end
end
"""
ldiv!(Y, A, B) -> Y
Compute `A \\ B` in-place and store the result in `Y`, returning the result.
The argument `A` should *not* be a matrix. Rather, instead of matrices it should be a
factorization object (e.g. produced by [`factorize`](@ref) or [`cholfact`](@ref)).
The reason for this is that factorization itself is both expensive and typically allocates memory
(although it can also be done in-place via, e.g., [`lufact!`](@ref)),
and performance-critical situations requiring `ldiv!` usually also require fine-grained
control over the factorization of `A`.
"""
ldiv!(Y, A, B)
"""
ldiv!(A, B)
Compute `A \\ B` in-place and overwriting `B` to store the result.
The argument `A` should *not* be a matrix. Rather, instead of matrices it should be a
factorization object (e.g. produced by [`factorize`](@ref) or [`cholfact`](@ref)).
The reason for this is that factorization itself is both expensive and typically allocates memory
(although it can also be done in-place via, e.g., [`lufact!`](@ref)),
and performance-critical situations requiring `ldiv!` usually also require fine-grained
control over the factorization of `A`.
"""
ldiv!(A, B)
"""
rdiv!(A, B)
Compute `A / B` in-place and overwriting `A` to store the result.
The argument `B` should *not* be a matrix. Rather, instead of matrices it should be a
factorization object (e.g. produced by [`factorize`](@ref) or [`cholfact`](@ref)).
The reason for this is that factorization itself is both expensive and typically allocates memory
(although it can also be done in-place via, e.g., [`lufact!`](@ref)),
and performance-critical situations requiring `rdiv!` usually also require fine-grained
control over the factorization of `B`.
"""
rdiv!(A, B)
copy_oftype(A::AbstractArray{T}, ::Type{T}) where {T} = copy(A)
copy_oftype(A::AbstractArray{T,N}, ::Type{S}) where {T,N,S} = convert(AbstractArray{S,N}, A)
include("adjtrans.jl")
include("transpose.jl")
include("conjarray.jl")
include("rowvector.jl")
include("exceptions.jl")
include("generic.jl")
include("blas.jl")
include("matmul.jl")
include("lapack.jl")
include("dense.jl")
include("tridiag.jl")
include("triangular.jl")
include("factorization.jl")
include("qr.jl")
include("hessenberg.jl")
include("lq.jl")
include("eigen.jl")
include("svd.jl")
include("symmetric.jl")
include("cholesky.jl")
include("lu.jl")
include("bunchkaufman.jl")
include("diagonal.jl")
include("bidiag.jl")
include("uniformscaling.jl")
include("givens.jl")
include("special.jl")
include("bitarray.jl")
include("ldlt.jl")
include("schur.jl")
include("deprecated.jl")
const ⋅ = dot
const × = cross
export ⋅, ×
function versioninfo(io::IO=STDOUT)
if Base.libblas_name == "libopenblas" || BLAS.vendor() == :openblas || BLAS.vendor() == :openblas64
openblas_config = BLAS.openblas_get_config()
println(io, "BLAS: libopenblas (", openblas_config, ")")
else
println(io, "BLAS: ",Base.libblas_name)
end
println(io, "LAPACK: ",Base.liblapack_name)
end
function __init__()
try
BLAS.check()
if BLAS.vendor() == :mkl
ccall((:MKL_Set_Interface_Layer, Base.libblas_name), Cvoid, (Cint,), USE_BLAS64 ? 1 : 0)
end
catch ex
Base.showerror_nostdio(ex,
"WARNING: Error during initialization of module LinearAlgebra")
end
end
end # module LinearAlgebra