This file is indexed.

/usr/share/julia/test/linalg/lu.jl is in julia-common 0.4.7-6.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
# This file is a part of Julia. License is MIT: http://julialang.org/license

debug = false

using Base.Test
import Base.LinAlg.BlasInt, Base.LinAlg.BlasFloat

n = 10

# Split n into 2 parts for tests needing two matrices
n1 = div(n, 2)
n2 = 2*n1

srand(1234321)

a = rand(n,n)

areal = randn(n,n)/2
aimg  = randn(n,n)/2
breal = randn(n,2)/2
bimg  = randn(n,2)/2
creal = randn(n)/2
cimg  = randn(n)/2
dureal = randn(n-1)/2
duimg  = randn(n-1)/2
dlreal = randn(n-1)/2
dlimg  = randn(n-1)/2
dreal = randn(n)/2
dimg  = randn(n)/2

for eltya in (Float32, Float64, Complex64, Complex128, BigFloat, Int)
    a = eltya == Int ? rand(1:7, n, n) : convert(Matrix{eltya}, eltya <: Complex ? complex(areal, aimg) : areal)
    d = eltya == Int ? Tridiagonal(rand(1:7, n-1), rand(1:7, n), rand(1:7, n-1)) : convert(Tridiagonal{eltya}, eltya <: Complex ? Tridiagonal(complex(dlreal, dlimg), complex(dreal, dimg), complex(dureal, duimg)) : Tridiagonal(dlreal, dreal, dureal))
    ε = εa = eps(abs(float(one(eltya))))

    if eltya <: BlasFloat
        num = rand(eltya)
        @test lu(num) == (one(eltya),num,1)
        @test_approx_eq full(lufact(num)) eltya[num]
    end
    for eltyb in (Float32, Float64, Complex64, Complex128, Int)
        b = eltyb == Int ? rand(1:5, n, 2) : convert(Matrix{eltyb}, eltyb <: Complex ? complex(breal, bimg) : breal)
        c = eltyb == Int ? rand(1:5, n) : convert(Vector{eltyb}, eltyb <: Complex ? complex(creal, cimg) : creal)
        εb = eps(abs(float(one(eltyb))))
        ε = max(εa,εb)

debug && println("(Automatic) Square LU decomposition")
        κ     = cond(a,1)
        lua   = factorize(a)
        @test_throws KeyError lua[:Z]
        l,u,p = lua[:L], lua[:U], lua[:p]
        ll,ul,pl = lu(a)
        @test_approx_eq ll * ul a[pl,:]
        @test_approx_eq l*u a[p,:]
        @test_approx_eq (l*u)[invperm(p),:] a
        @test_approx_eq a * inv(lua) eye(n)
        @test norm(a*(lua\b) - b, 1) < ε*κ*n*2 # Two because the right hand side has two columns
        @test norm(a'*(lua'\b) - b, 1) < ε*κ*n*2 # Two because the right hand side has two columns
        @test norm(a'*(lua'\a') - a', 1) < ε*κ*n^2
        @test norm(a*(lua\c) - c, 1) < ε*κ*n # c is a vector
        @test norm(a'*(lua'\c) - c, 1) < ε*κ*n # c is a vector
        @test_approx_eq full(lua) a
        if eltya <: Real && eltyb <: Real
            @test norm(a.'*(lua.'\b) - b,1) < ε*κ*n*2 # Two because the right hand side has two columns
            @test norm(a.'*(lua.'\c) - c,1) < ε*κ*n
        end
        if eltya <: BlasFloat && eltyb <: BlasFloat
            e = rand(eltyb,n,n)
            @test norm(e/lua - e/a,1) < ε*κ*n^2
        end

debug && println("Tridiagonal LU")
        κd    = cond(full(d),1)
        lud   = lufact(d)
        @test lufact(lud) == lud
        @test_throws KeyError lud[:Z]
        @test_approx_eq lud[:L]*lud[:U] lud[:P]*full(d)
        @test_approx_eq lud[:L]*lud[:U] full(d)[lud[:p],:]
        @test_approx_eq full(lud) d
        @test norm(d*(lud\b) - b, 1) < ε*κd*n*2 # Two because the right hand side has two columns
        if eltya <: Real
            @test norm((lud.'\b) - full(d.')\b, 1) < ε*κd*n*2 # Two because the right hand side has two columns
        end
        if eltya <: Complex
            @test norm((lud'\b) - full(d')\b, 1) < ε*κd*n*2 # Two because the right hand side has two columns
        end
        f = zeros(eltyb, n+1)
        @test_throws DimensionMismatch lud\f
        @test_throws DimensionMismatch lud.'\f
        @test_throws DimensionMismatch lud'\f

        if eltya <: BlasFloat && eltyb <: BlasFloat
            e = rand(eltyb,n,n)
            @test norm(e/lud - e/d,1) < ε*κ*n^2
            @test norm((lud.'\e') - full(d.')\e',1) < ε*κd*n^2
            #test singular
            du = rand(eltya,n-1)
            dl = rand(eltya,n-1)
            dd = rand(eltya,n)
            dd[1] = zero(eltya)
            du[1] = zero(eltya)
            dl[1] = zero(eltya)
            zT = Tridiagonal(dl,dd,du)
            @test lufact(zT).info == 1
        end

debug && println("Thin LU")
        lua   = @inferred lufact(a[:,1:n1])
        @test_approx_eq lua[:L]*lua[:U] lua[:P]*a[:,1:n1]

debug && println("Fat LU")
        lua   = lufact(a[1:n1,:])
        @test_approx_eq lua[:L]*lua[:U] lua[:P]*a[1:n1,:]
    end
end

# test conversion routine
a = Tridiagonal(rand(9),rand(10),rand(9))
fa = full(a)
falu = lufact(fa)
alu = lufact(a)
falu = convert(typeof(falu),alu)
@test full(alu) == fa

# Test rational matrices
## Integrate in general tests when more linear algebra is implemented in julia
a = convert(Matrix{Rational{BigInt}}, rand(1:10//1,n,n))/n
b = rand(1:10,n,2)
@inferred lufact(a)
lua   = factorize(a)
l,u,p = lua[:L], lua[:U], lua[:p]
@test_approx_eq l*u a[p,:]
@test_approx_eq l[invperm(p),:]*u a
@test_approx_eq a * inv(lua) eye(n)
@test_approx_eq a*(lua\b) b
@test_approx_eq @inferred(det(a)) det(Array{Float64}(a))
## Hilbert Matrix (very ill conditioned)
## Testing Rational{BigInt} and BigFloat version
nHilbert = 50
H = Rational{BigInt}[1//(i+j-1) for i = 1:nHilbert,j = 1:nHilbert]
Hinv = Rational{BigInt}[(-1)^(i+j)*(i+j-1)*binomial(nHilbert+i-1,nHilbert-j)*binomial(nHilbert+j-1,nHilbert-i)*binomial(i+j-2,i-1)^2 for i = big(1):nHilbert,j=big(1):nHilbert]
@test inv(H) == Hinv
with_bigfloat_precision(2^10) do
    @test norm(Array{Float64}(inv(float(H)) - float(Hinv))) < 1e-100
end

# Test balancing in eigenvector calculations
for elty in (Float32, Float64, Complex64, Complex128)
    A = convert(Matrix{elty}, [ 3.0     -2.0      -0.9     2*eps(real(one(elty)));
                               -2.0      4.0       1.0    -eps(real(one(elty)));
                               -eps(real(one(elty)))/4  eps(real(one(elty)))/2  -1.0     0;
                               -0.5     -0.5       0.1     1.0])
    F = eigfact(A, permute=false, scale=false)
    eig(A, permute=false, scale=false)
    @test_approx_eq F[:vectors]*Diagonal(F[:values])/F[:vectors] A
    F = eigfact(A)
    # @test norm(F[:vectors]*Diagonal(F[:values])/F[:vectors] - A) > 0.01
end

@test @inferred(logdet(Complex64[1.0f0 0.5f0; 0.5f0 -1.0f0])) === 0.22314355f0 + 3.1415927f0im
@test_throws DomainError logdet([1 1; 1 -1])