-
-
Notifications
You must be signed in to change notification settings - Fork 5.5k
/
Copy pathdense.jl
467 lines (416 loc) · 13.8 KB
/
dense.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
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
# Linear algebra functions for dense matrices in column major format
## BLAS cutoff threshold constants
const SCAL_CUTOFF = 2048
const DOT_CUTOFF = 128
const ASUM_CUTOFF = 32
const NRM2_CUTOFF = 32
function scale!{T<:BlasFloat}(X::Array{T}, s::T)
if length(X) < SCAL_CUTOFF
generic_scale!(X, s)
else
BLAS.scal!(length(X), s, X, 1)
end
X
end
scale!{T<:BlasFloat}(X::Array{T}, s::Number) = scale!(X, convert(T, s))
function scale!{T<:BlasComplex}(X::Array{T}, s::Real)
R = typeof(real(zero(T)))
BLAS.scal!(2*length(X), convert(R,s), convert(Ptr{R},pointer(X)), 1)
X
end
#Test whether a matrix is positive-definite
isposdef!{T<:BlasFloat}(A::StridedMatrix{T}, UL::Symbol) = LAPACK.potrf!(char_uplo(UL), A)[2] == 0
isposdef!(A::StridedMatrix) = ishermitian(A) && isposdef!(A, :U)
isposdef{T}(A::AbstractMatrix{T}, UL::Symbol) = (S = typeof(sqrt(one(T))); isposdef!(S == T ? copy(A) : convert(AbstractMatrix{S}, A), UL))
isposdef{T}(A::AbstractMatrix{T}) = (S = typeof(sqrt(one(T))); isposdef!(S == T ? copy(A) : convert(AbstractMatrix{S}, A)))
isposdef(x::Number) = imag(x)==0 && real(x) > 0
stride1(x::Array) = 1
stride1(x::StridedVector) = stride(x, 1)::Int
import Base: mapreduce_seq_impl, AbsFun, Abs2Fun, AddFun
mapreduce_seq_impl{T<:BlasReal}(::AbsFun, ::AddFun, a::Union(Array{T},StridedVector{T}), ifirst::Int, ilast::Int) =
BLAS.asum(ilast-ifirst+1, pointer(a, ifirst), stride1(a))
function mapreduce_seq_impl{T<:BlasReal}(::Abs2Fun, ::AddFun, a::Union(Array{T},StridedVector{T}), ifirst::Int, ilast::Int)
n = ilast-ifirst+1
px = pointer(a, ifirst)
incx = stride1(a)
BLAS.dot(n, px, incx, px, incx)
end
function mapreduce_seq_impl{T<:BlasComplex}(::Abs2Fun, ::AddFun, a::Union(Array{T},StridedVector{T}), ifirst::Int, ilast::Int)
n = ilast-ifirst+1
px = pointer(a, ifirst)
incx = stride1(a)
real(BLAS.dotc(n, px, incx, px, incx))
end
function norm{T<:BlasFloat, TI<:Integer}(x::StridedVector{T}, rx::Union(UnitRange{TI},Range{TI}))
(minimum(rx) < 1 || maximum(rx) > length(x)) && throw(BoundsError())
BLAS.nrm2(length(rx), pointer(x)+(first(rx)-1)*sizeof(T), step(rx))
end
vecnorm1{T<:BlasReal}(x::Union(Array{T},StridedVector{T})) =
length(x) < ASUM_CUTOFF ? generic_vecnorm1(x) : BLAS.asum(x)
vecnorm2{T<:BlasFloat}(x::Union(Array{T},StridedVector{T})) =
length(x) < NRM2_CUTOFF ? generic_vecnorm2(x) : BLAS.nrm2(x)
function triu!(M::AbstractMatrix, k::Integer)
m, n = size(M)
idx = 1
for j = 0:n-1
ii = min(max(0, j+1-k), m)
for i = (idx+ii):(idx+m-1)
M[i] = zero(M[i])
end
idx += m
end
M
end
triu(M::Matrix, k::Integer) = triu!(copy(M), k)
function tril!(M::AbstractMatrix, k::Integer)
m, n = size(M)
idx = 1
for j = 0:n-1
ii = min(max(0, j-k), m)
for i = idx:(idx+ii-1)
M[i] = zero(M[i])
end
idx += m
end
M
end
tril(M::Matrix, k::Integer) = tril!(copy(M), k)
function gradient(F::Vector, h::Vector)
n = length(F)
T = typeof(one(eltype(F))/one(eltype(h)))
g = Array(T,n)
if n == 1
g[1] = zero(T)
elseif n > 1
g[1] = (F[2] - F[1]) / (h[2] - h[1])
g[n] = (F[n] - F[n-1]) / (h[end] - h[end-1])
if n > 2
h = h[3:n] - h[1:n-2]
g[2:n-1] = (F[3:n] - F[1:n-2]) ./ h
end
end
g
end
function diagind(m::Integer, n::Integer, k::Integer=0)
-m <= k <= n || throw(BoundsError())
k <= 0 ? range(1-k, m+1, min(m+k, n)) : range(k*m+1, m+1, min(m, n-k))
end
diagind(A::AbstractMatrix, k::Integer=0) = diagind(size(A,1), size(A,2), k)
diag(A::AbstractMatrix, k::Integer=0) = A[diagind(A,k)]
function diagm{T}(v::AbstractVector{T}, k::Integer=0)
n = length(v) + abs(k)
A = zeros(T,n,n)
A[diagind(A,k)] = v
A
end
diagm(x::Number) = (X = Array(typeof(x),1,1); X[1,1] = x; X)
function trace{T}(A::Matrix{T})
n = chksquare(A)
t = zero(T)
for i=1:n
t += A[i,i]
end
t
end
function kron{T,S}(a::Matrix{T}, b::Matrix{S})
R = Array(promote_type(T,S), size(a,1)*size(b,1), size(a,2)*size(b,2))
m = 1
for j = 1:size(a,2), l = 1:size(b,2), i = 1:size(a,1)
aij = a[i,j]
for k = 1:size(b,1)
R[m] = aij*b[k,l]
m += 1
end
end
R
end
kron(a::Number, b::Union(Number, Vector, Matrix)) = a * b
kron(a::Union(Vector, Matrix), b::Number) = a * b
kron(a::Vector, b::Vector)=vec(kron(reshape(a,length(a),1),reshape(b,length(b),1)))
kron(a::Matrix, b::Vector)=kron(a,reshape(b,length(b),1))
kron(a::Vector, b::Matrix)=kron(reshape(a,length(a),1),b)
^(A::Matrix, p::Integer) = p < 0 ? inv(A^-p) : Base.power_by_squaring(A,p)
function ^(A::Matrix, p::Number)
isinteger(p) && return A^integer(real(p))
chksquare(A)
v, X = eig(A)
any(v.<0) && (v = complex(v))
Xinv = ishermitian(A) ? X' : inv(X)
scale(X, v.^p)*Xinv
end
# Matrix exponential
expm{T<:BlasFloat}(A::StridedMatrix{T}) = expm!(copy(A))
expm{T<:Integer}(A::StridedMatrix{T}) = expm!(float(A))
expm(x::Number) = exp(x)
## Destructive matrix exponential using algorithm from Higham, 2008,
## "Functions of Matrices: Theory and Computation", SIAM
function expm!{T<:BlasFloat}(A::StridedMatrix{T})
n = chksquare(A)
n<2 && return exp(A)
ilo, ihi, scale = LAPACK.gebal!('B', A) # modifies A
nA = norm(A, 1)
I = eye(T,n)
## For sufficiently small nA, use lower order Padé-Approximations
if (nA <= 2.1)
if nA > 0.95
C = T[17643225600.,8821612800.,2075673600.,302702400.,
30270240., 2162160., 110880., 3960.,
90., 1.]
elseif nA > 0.25
C = T[17297280.,8648640.,1995840.,277200.,
25200., 1512., 56., 1.]
elseif nA > 0.015
C = T[30240.,15120.,3360.,
420., 30., 1.]
else
C = T[120.,60.,12.,1.]
end
A2 = A * A
P = copy(I)
U = C[2] * P
V = C[1] * P
for k in 1:(div(size(C, 1), 2) - 1)
k2 = 2 * k
P *= A2
U += C[k2 + 2] * P
V += C[k2 + 1] * P
end
U = A * U
X = V + U
LAPACK.gesv!(V-U, X)
else
s = log2(nA/5.4) # power of 2 later reversed by squaring
if s > 0
si = ceil(Int,s)
A /= convert(T,2^si)
end
CC = T[64764752532480000.,32382376266240000.,7771770303897600.,
1187353796428800., 129060195264000., 10559470521600.,
670442572800., 33522128640., 1323241920.,
40840800., 960960., 16380.,
182., 1.]
A2 = A * A
A4 = A2 * A2
A6 = A2 * A4
U = A * (A6 * (CC[14]*A6 + CC[12]*A4 + CC[10]*A2) +
CC[8]*A6 + CC[6]*A4 + CC[4]*A2 + CC[2]*I)
V = A6 * (CC[13]*A6 + CC[11]*A4 + CC[9]*A2) +
CC[7]*A6 + CC[5]*A4 + CC[3]*A2 + CC[1]*I
X = V + U
LAPACK.gesv!(V-U, X)
if s > 0 # squaring to reverse dividing by power of 2
for t=1:si X *= X end
end
end
# Undo the balancing
for j = ilo:ihi
scj = scale[j]
for i = 1:n
X[j,i] *= scj
end
for i = 1:n
X[i,j] /= scj
end
end
if ilo > 1 # apply lower permutations in reverse order
for j in (ilo-1):-1:1 rcswap!(j, int(scale[j]), X) end
end
if ihi < n # apply upper permutations in forward order
for j in (ihi+1):n rcswap!(j, int(scale[j]), X) end
end
X
end
## Swap rows i and j and columns i and j in X
function rcswap!{T<:Number}(i::Integer, j::Integer, X::StridedMatrix{T})
for k = 1:size(X,1)
X[k,i], X[k,j] = X[k,j], X[k,i]
end
for k = 1:size(X,2)
X[i,k], X[j,k] = X[j,k], X[i,k]
end
end
function sqrtm{T<:Real}(A::StridedMatrix{T})
issym(A) && return sqrtm(Symmetric(A))
n = chksquare(A)
SchurF = schurfact(complex(A))
R = full(sqrtm(UpperTriangular(SchurF[:T])))
retmat = SchurF[:vectors]*R*SchurF[:vectors]'
all(imag(retmat) .== 0) ? real(retmat) : retmat
end
function sqrtm{T<:Complex}(A::StridedMatrix{T})
ishermitian(A) && return sqrtm(Hermitian(A))
n = chksquare(A)
SchurF = schurfact(A)
R = full(sqrtm(UpperTriangular(SchurF[:T])))
SchurF[:vectors]*R*SchurF[:vectors]'
end
sqrtm(a::Number) = (b = sqrt(complex(a)); imag(b) == 0 ? real(b) : b)
sqrtm(a::Complex) = sqrt(a)
function inv{S}(A::StridedMatrix{S})
T = typeof(one(S)/one(S))
Ac = convert(AbstractMatrix{T}, A)
if istriu(Ac)
Ai = inv(UpperTriangular(A))
elseif istril(Ac)
Ai = inv(LowerTriangular(A))
else
Ai = inv(lufact(Ac))
end
return convert(typeof(Ac), Ai)
end
function factorize{T}(A::Matrix{T})
m, n = size(A)
if m == n
if m == 1 return A[1] end
utri = true
utri1 = true
herm = true
sym = true
for j = 1:n-1, i = j+1:m
if utri1
if A[i,j] != 0
utri1 = i == j + 1
utri = false
end
end
if sym
sym &= A[i,j] == A[j,i]
end
if herm
herm &= A[i,j] == conj(A[j,i])
end
if !(utri1|herm|sym) break end
end
ltri = true
ltri1 = true
for j = 3:n, i = 1:j-2
ltri1 &= A[i,j] == 0
if !ltri1 break end
end
if ltri1
for i = 1:n-1
if A[i,i+1] != 0
ltri &= false
break
end
end
if ltri
if utri
return Diagonal(A)
end
if utri1
return Bidiagonal(diag(A), diag(A, -1), false)
end
return LowerTriangular(A)
end
if utri
return Bidiagonal(diag(A), diag(A, 1), true)
end
if utri1
if (herm & (T <: Complex)) | sym
try
return ldltfact!(SymTridiagonal(diag(A), diag(A, -1)))
end
end
return lufact(Tridiagonal(diag(A, -1), diag(A), diag(A, 1)))
end
end
if utri
return UpperTriangular(A)
end
if herm
try
return cholfact(A)
end
return factorize(Hermitian(A))
end
if sym
return factorize(Symmetric(A))
end
return lufact(A)
end
qrfact(A,typeof(zero(T)/sqrt(zero(T) + zero(T)))<:BlasFloat?Val{true}:Val{false}) # Generic pivoted QR not implemented yet
end
(\)(a::Vector, B::StridedVecOrMat) = (\)(reshape(a, length(a), 1), B)
for (T1,PIVOT) in ((BlasFloat,Val{true}),(Any,Val{false}))
@eval function (\){T<:$T1}(A::StridedMatrix{T}, B::StridedVecOrMat)
m, n = size(A)
if m == n
if istril(A)
return istriu(A) ? \(Diagonal(A),B) : \(LowerTriangular(A),B)
end
istriu(A) && return \(UpperTriangular(A),B)
return \(lufact(A),B)
end
return qrfact(A,$PIVOT)\B
end
end
## Moore-Penrose pseudoinverse
function pinv{T}(A::StridedMatrix{T}, tol::Real)
m, n = size(A)
(m == 0 || n == 0) && return Array(T, n, m)
if istril(A)
if( istriu(A) )
maxabsA = maximum(abs(diag(A)))
B = zeros(T,n,m);
for i = 1:min(m,n)
if( abs(A[i,i]) > tol*maxabsA && isfinite(one(T)/A[i,i]) )
B[i,i] = one(T)/A[i,i]
end
end
return B;
end
end
SVD = svdfact(A, thin=true)
S = eltype(SVD[:S])
Sinv = zeros(S, length(SVD[:S]))
index = SVD[:S] .> tol*maximum(SVD[:S])
Sinv[index] = one(S) ./ SVD[:S][index]
Sinv[find(!isfinite(Sinv))] = zero(S)
return SVD[:Vt]'scale(Sinv, SVD[:U]')
end
function pinv{T}(A::StridedMatrix{T})
tol = eps(real(float(one(T))))*maximum(size(A))
return pinv(A, tol)
end
pinv(a::StridedVector) = pinv(reshape(a, length(a), 1))
pinv(x::Number) = isfinite(one(x)/x) ? one(x)/x : zero(x)
## Basis for null space
function nullspace{T}(A::StridedMatrix{T})
m, n = size(A)
(m == 0 || n == 0) && return eye(T, n)
SVD = svdfact(A, thin=false)
indstart = sum(SVD[:S] .> max(m,n)*maximum(SVD[:S])*eps(eltype(SVD[:S]))) + 1
return SVD[:V][:,indstart:end]
end
nullspace(a::StridedVector) = nullspace(reshape(a, length(a), 1))
function cond(A::AbstractMatrix, p::Real=2)
if p == 2
v = svdvals(A)
maxv = maximum(v)
return maxv == 0.0 ? oftype(real(A[1,1]),Inf) : maxv / minimum(v)
elseif p == 1 || p == Inf
chksquare(A)
return cond(lufact(A), p)
end
throw(ArgumentError("p-norm must be 1, 2 or Inf, got $p"))
end
## Lyapunov and Sylvester equation
# AX + XB + C = 0
function sylvester{T<:BlasFloat}(A::StridedMatrix{T},B::StridedMatrix{T},C::StridedMatrix{T})
RA, QA = schur(A)
RB, QB = schur(B)
D = -Ac_mul_B(QA,C*QB)
Y, scale = LAPACK.trsyl!('N','N', RA, RB, D)
scale!(QA*A_mul_Bc(Y,QB), inv(scale))
end
sylvester{T<:Integer}(A::StridedMatrix{T},B::StridedMatrix{T},C::StridedMatrix{T}) = sylvester(float(A), float(B), float(C))
# AX + XA' + C = 0
function lyap{T<:BlasFloat}(A::StridedMatrix{T},C::StridedMatrix{T})
R, Q = schur(A)
D = -Ac_mul_B(Q,C*Q)
Y, scale = LAPACK.trsyl!('N', T <: Complex ? 'C' : 'T', R, R, D)
scale!(Q*A_mul_Bc(Y,Q), inv(scale))
end
lyap{T<:Integer}(A::StridedMatrix{T},C::StridedMatrix{T}) = lyap(float(A), float(C))
lyap{T<:Number}(a::T, c::T) = -c/(2a)