• Home
  • Features
  • Pricing
  • Docs
  • Announcements
  • Sign In

JuliaLang / julia / #37933

15 Oct 2024 07:30AM UTC coverage: 87.724% (+1.3%) from 86.441%
#37933

push

local

web-flow
Fix markdown list in installation.md (#56165)

Documenter.jl requires all trailing list content to follow the same
indentation as the header. So, in the current view
(https://docs.julialang.org/en/v1/manual/installation/#Command-line-arguments)
the list appears broken.

78955 of 90004 relevant lines covered (87.72%)

16956008.11 hits per line

Source File
Press 'n' to go to next uncovered line, 'b' for previous

93.05
/stdlib/LinearAlgebra/src/triangular.jl
1
# This file is a part of Julia. License is MIT: https://julialang.org/license
2

3
## Triangular
4

5
# could be renamed to Triangular when that name has been fully deprecated
6
"""
7
    AbstractTriangular
8

9
Supertype of triangular matrix types such as [`LowerTriangular`](@ref), [`UpperTriangular`](@ref),
10
[`UnitLowerTriangular`](@ref) and [`UnitUpperTriangular`](@ref).
11
"""
12
abstract type AbstractTriangular{T} <: AbstractMatrix{T} end
13

14
# First loop through all methods that don't need special care for upper/lower and unit diagonal
15
for t in (:LowerTriangular, :UnitLowerTriangular, :UpperTriangular, :UnitUpperTriangular)
16
    @eval begin
17
        struct $t{T,S<:AbstractMatrix{T}} <: AbstractTriangular{T}
18
            data::S
19

20
            function $t{T,S}(data) where {T,S<:AbstractMatrix{T}}
16✔
21
                require_one_based_indexing(data)
96,959✔
22
                checksquare(data)
96,960✔
23
                new{T,S}(data)
96,958✔
24
            end
25
        end
26
        $t(A::$t) = A
11,958✔
27
        $t{T}(A::$t{T}) where {T} = A
60✔
28
        $t(A::AbstractMatrix) = $t{eltype(A), typeof(A)}(A)
96,925✔
29
        $t{T}(A::AbstractMatrix) where {T} = $t(convert(AbstractMatrix{T}, A))
28✔
30
        $t{T}(A::$t) where {T} = $t(convert(AbstractMatrix{T}, A.data))
3,624✔
31

32
        AbstractMatrix{T}(A::$t) where {T} = $t{T}(A)
3,592✔
33
        AbstractMatrix{T}(A::$t{T}) where {T} = copy(A)
×
34

35
        size(A::$t) = size(A.data)
480,864✔
36
        axes(A::$t) = axes(A.data)
780,182✔
37

38
        # For A<:AbstractTriangular, similar(A[, neweltype]) should yield a matrix with the same
39
        # triangular type and underlying storage type as A. The following method covers these cases.
40
        similar(A::$t, ::Type{T}) where {T} = $t(similar(parent(A), T))
5,777✔
41
        # On the other hand, similar(A, [neweltype,] shape...) should yield a matrix of the underlying
42
        # storage type of A (not wrapped in a triangular type). The following method covers these cases.
43
        similar(A::$t, ::Type{T}, dims::Dims{N}) where {T,N} = similar(parent(A), T, dims)
32,973✔
44

45
        copy(A::$t) = $t(copy(A.data))
1✔
46
        Base.unaliascopy(A::$t) = $t(Base.unaliascopy(A.data))
2,284✔
47

48
        real(A::$t{<:Complex}) = (B = real(A.data); $t(B))
×
49
        real(A::$t{<:Complex, <:StridedMaybeAdjOrTransMat}) = $t(real.(A))
370✔
50
    end
51
end
52

53
"""
54
    LowerTriangular(A::AbstractMatrix)
55

56
Construct a `LowerTriangular` view of the matrix `A`.
57

58
# Examples
59
```jldoctest
60
julia> A = [1.0 2.0 3.0; 4.0 5.0 6.0; 7.0 8.0 9.0]
61
3×3 Matrix{Float64}:
62
 1.0  2.0  3.0
63
 4.0  5.0  6.0
64
 7.0  8.0  9.0
65

66
julia> LowerTriangular(A)
67
3×3 LowerTriangular{Float64, Matrix{Float64}}:
68
 1.0   ⋅    ⋅
69
 4.0  5.0   ⋅
70
 7.0  8.0  9.0
71
```
72
"""
73
LowerTriangular
74
"""
75
    UpperTriangular(A::AbstractMatrix)
76

77
Construct an `UpperTriangular` view of the matrix `A`.
78

79
# Examples
80
```jldoctest
81
julia> A = [1.0 2.0 3.0; 4.0 5.0 6.0; 7.0 8.0 9.0]
82
3×3 Matrix{Float64}:
83
 1.0  2.0  3.0
84
 4.0  5.0  6.0
85
 7.0  8.0  9.0
86

87
julia> UpperTriangular(A)
88
3×3 UpperTriangular{Float64, Matrix{Float64}}:
89
 1.0  2.0  3.0
90
  ⋅   5.0  6.0
91
  ⋅    ⋅   9.0
92
```
93
"""
94
UpperTriangular
95
"""
96
    UnitLowerTriangular(A::AbstractMatrix)
97

98
Construct a `UnitLowerTriangular` view of the matrix `A`.
99
Such a view has the [`oneunit`](@ref) of the [`eltype`](@ref)
100
of `A` on its diagonal.
101

102
# Examples
103
```jldoctest
104
julia> A = [1.0 2.0 3.0; 4.0 5.0 6.0; 7.0 8.0 9.0]
105
3×3 Matrix{Float64}:
106
 1.0  2.0  3.0
107
 4.0  5.0  6.0
108
 7.0  8.0  9.0
109

110
julia> UnitLowerTriangular(A)
111
3×3 UnitLowerTriangular{Float64, Matrix{Float64}}:
112
 1.0   ⋅    ⋅
113
 4.0  1.0   ⋅
114
 7.0  8.0  1.0
115
```
116
"""
117
UnitLowerTriangular
118
"""
119
    UnitUpperTriangular(A::AbstractMatrix)
120

121
Construct an `UnitUpperTriangular` view of the matrix `A`.
122
Such a view has the [`oneunit`](@ref) of the [`eltype`](@ref)
123
of `A` on its diagonal.
124

125
# Examples
126
```jldoctest
127
julia> A = [1.0 2.0 3.0; 4.0 5.0 6.0; 7.0 8.0 9.0]
128
3×3 Matrix{Float64}:
129
 1.0  2.0  3.0
130
 4.0  5.0  6.0
131
 7.0  8.0  9.0
132

133
julia> UnitUpperTriangular(A)
134
3×3 UnitUpperTriangular{Float64, Matrix{Float64}}:
135
 1.0  2.0  3.0
136
  ⋅   1.0  6.0
137
  ⋅    ⋅   1.0
138
```
139
"""
140
UnitUpperTriangular
141

142
const UpperOrUnitUpperTriangular{T,S} = Union{UpperTriangular{T,S}, UnitUpperTriangular{T,S}}
143
const LowerOrUnitLowerTriangular{T,S} = Union{LowerTriangular{T,S}, UnitLowerTriangular{T,S}}
144
const UpperOrLowerTriangular{T,S} = Union{UpperOrUnitUpperTriangular{T,S}, LowerOrUnitLowerTriangular{T,S}}
145

146
uppertriangular(M) = UpperTriangular(M)
1,089✔
147
lowertriangular(M) = LowerTriangular(M)
262✔
148

149
uppertriangular(U::UpperOrUnitUpperTriangular) = U
17✔
150
lowertriangular(U::LowerOrUnitLowerTriangular) = U
1✔
151

152
Base.dataids(A::UpperOrLowerTriangular) = Base.dataids(A.data)
24,970✔
153

154
imag(A::UpperTriangular) = UpperTriangular(imag(A.data))
1✔
155
imag(A::LowerTriangular) = LowerTriangular(imag(A.data))
1✔
156
imag(A::UpperTriangular{<:Any,<:StridedMaybeAdjOrTransMat}) = imag.(A)
12✔
157
imag(A::LowerTriangular{<:Any,<:StridedMaybeAdjOrTransMat}) = imag.(A)
10✔
158
function imag(A::UnitLowerTriangular)
9✔
159
    L = LowerTriangular(A.data)
9✔
160
    Lim = similar(L) # must be mutable to set diagonals to zero
9✔
161
    Lim .= imag.(L)
9✔
162
    for i in 1:size(Lim,1)
9✔
163
        Lim[i,i] = zero(Lim[i,i])
67✔
164
    end
125✔
165
    return Lim
9✔
166
end
167
function imag(A::UnitUpperTriangular)
9✔
168
    U = UpperTriangular(A.data)
9✔
169
    Uim = similar(U) # must be mutable to set diagonals to zero
9✔
170
    Uim .= imag.(U)
9✔
171
    for i in 1:size(Uim,1)
9✔
172
        Uim[i,i] = zero(Uim[i,i])
67✔
173
    end
125✔
174
    return Uim
9✔
175
end
176

177
parent(A::UpperOrLowerTriangular) = A.data
179,645✔
178

179
# For strided matrices, we may only loop over the filled triangle
180
copy(A::UpperOrLowerTriangular{<:Any, <:StridedMaybeAdjOrTransMat}) = copyto!(similar(A), A)
4,199✔
181

182
# then handle all methods that requires specific handling of upper/lower and unit diagonal
183

184
function full!(A::LowerTriangular)
7✔
185
    B = A.data
7✔
186
    tril!(B)
7✔
187
    B
7✔
188
end
189
function full!(A::UnitLowerTriangular)
7✔
190
    B = A.data
7✔
191
    tril!(B)
7✔
192
    for i = 1:size(A,1)
7✔
193
        B[i,i] = oneunit(eltype(B))
63✔
194
    end
119✔
195
    B
7✔
196
end
197
function full!(A::UpperTriangular)
7✔
198
    B = A.data
7✔
199
    triu!(B)
7✔
200
    B
7✔
201
end
202
function full!(A::UnitUpperTriangular)
7✔
203
    B = A.data
7✔
204
    triu!(B)
7✔
205
    for i = 1:size(A,1)
7✔
206
        B[i,i] = oneunit(eltype(B))
63✔
207
    end
119✔
208
    B
7✔
209
end
210

211
_shouldforwardindex(U::UpperTriangular, row::Integer, col::Integer) = row <= col
3,745,050✔
212
_shouldforwardindex(U::LowerTriangular, row::Integer, col::Integer) = row >= col
3,211,682✔
213
_shouldforwardindex(U::UnitUpperTriangular, row::Integer, col::Integer) = row < col
3,372,672✔
214
_shouldforwardindex(U::UnitLowerTriangular, row::Integer, col::Integer) = row > col
3,560,783✔
215

216
Base.isassigned(A::UpperOrLowerTriangular, i::Int, j::Int) =
630,900✔
217
    _shouldforwardindex(A, i, j) ? isassigned(A.data, i, j) : true
218

219
Base.isstored(A::UpperOrLowerTriangular, i::Int, j::Int) =
×
220
    _shouldforwardindex(A, i, j) ? Base.isstored(A.data, i, j) : false
221

222
@propagate_inbounds getindex(A::Union{UnitLowerTriangular{T}, UnitUpperTriangular{T}}, i::Int, j::Int) where {T} =
6,396,713✔
223
    _shouldforwardindex(A, i, j) ? A.data[i,j] : ifelse(i == j, oneunit(T), zero(T))
224
@propagate_inbounds getindex(A::Union{LowerTriangular, UpperTriangular}, i::Int, j::Int) =
6,959,597✔
225
    _shouldforwardindex(A, i, j) ? A.data[i,j] : _zero(A.data,j,i)
226

227
_shouldforwardindex(U::UpperTriangular, b::BandIndex) = b.band >= 0
5✔
228
_shouldforwardindex(U::LowerTriangular, b::BandIndex) = b.band <= 0
5✔
229
_shouldforwardindex(U::UnitUpperTriangular, b::BandIndex) = b.band > 0
8✔
230
_shouldforwardindex(U::UnitLowerTriangular, b::BandIndex) = b.band < 0
8✔
231

232
# these specialized getindex methods enable constant-propagation of the band
233
Base.@constprop :aggressive @propagate_inbounds function getindex(A::Union{UnitLowerTriangular{T}, UnitUpperTriangular{T}}, b::BandIndex) where {T}
234
    _shouldforwardindex(A, b) ? A.data[b] : ifelse(b.band == 0, oneunit(T), zero(T))
16✔
235
end
236
Base.@constprop :aggressive @propagate_inbounds function getindex(A::Union{LowerTriangular, UpperTriangular}, b::BandIndex)
237
    _shouldforwardindex(A, b) ? A.data[b] : _zero(A.data, b)
10✔
238
end
239

240
_zero_triangular_half_str(::Type{<:UpperOrUnitUpperTriangular}) = "lower"
510✔
241
_zero_triangular_half_str(::Type{<:LowerOrUnitLowerTriangular}) = "upper"
509✔
242

243
@noinline function throw_nonzeroerror(T, @nospecialize(x), i, j)
1,019✔
244
    Ts = _zero_triangular_half_str(T)
1,019✔
245
    Tn = nameof(T)
1,019✔
246
    throw(ArgumentError(
1,019✔
247
        lazy"cannot set index in the $Ts triangular part ($i, $j) of an $Tn matrix to a nonzero value ($x)"))
248
end
249
@noinline function throw_nononeerror(T, @nospecialize(x), i, j)
128✔
250
    Tn = nameof(T)
256✔
251
    throw(ArgumentError(
128✔
252
        lazy"cannot set index on the diagonal ($i, $j) of an $Tn matrix to a non-unit value ($x)"))
253
end
254

255
@propagate_inbounds function setindex!(A::UpperTriangular, x, i::Integer, j::Integer)
820✔
256
    if i > j
30,948✔
257
        iszero(x) || throw_nonzeroerror(typeof(A), x, i, j)
3,362✔
258
    else
259
        A.data[i,j] = x
27,843✔
260
    end
261
    return A
30,691✔
262
end
263

264
@propagate_inbounds function setindex!(A::UnitUpperTriangular, x, i::Integer, j::Integer)
882✔
265
    if i > j
1,658✔
266
        iszero(x) || throw_nonzeroerror(typeof(A), x, i, j)
1,101✔
267
    elseif i == j
810✔
268
        x == oneunit(x) || throw_nononeerror(typeof(A), x, i, j)
216✔
269
    else
270
        A.data[i,j] = x
594✔
271
    end
272
    return A
1,342✔
273
end
274

275
@propagate_inbounds function setindex!(A::LowerTriangular, x, i::Integer, j::Integer)
820✔
276
    if i < j
3,660✔
277
        iszero(x) || throw_nonzeroerror(typeof(A), x, i, j)
1,174✔
278
    else
279
        A.data[i,j] = x
2,742✔
280
    end
281
    return A
3,404✔
282
end
283

284
@propagate_inbounds function setindex!(A::UnitLowerTriangular, x, i::Integer, j::Integer)
882✔
285
    if i < j
1,659✔
286
        iszero(x) || throw_nonzeroerror(typeof(A), x, i, j)
1,100✔
287
    elseif i == j
812✔
288
        x == oneunit(x) || throw_nononeerror(typeof(A), x, i, j)
216✔
289
    else
290
        A.data[i,j] = x
596✔
291
    end
292
    return A
1,343✔
293
end
294

295
@noinline function throw_setindex_structuralzero_error(T, @nospecialize(x))
×
296
    Ts = _zero_triangular_half_str(T)
×
297
    Tn = nameof(T)
×
298
    throw(ArgumentError(
×
299
        lazy"cannot set indices in the $Ts triangular part of an $Tn matrix to a nonzero value ($x)"))
×
300
end
301

302
@inline function fill!(A::UpperTriangular, x)
303
    iszero(x) || throw_setindex_structuralzero_error(typeof(A), x)
109✔
304
    for col in axes(A,2), row in firstindex(A,1):col
109✔
305
        @inbounds A.data[row, col] = x
931✔
306
    end
1,192✔
307
    A
109✔
308
end
309
@inline function fill!(A::LowerTriangular, x)
310
    iszero(x) || throw_setindex_structuralzero_error(typeof(A), x)
114✔
311
    for col in axes(A,2), row in col:lastindex(A,1)
114✔
312
        @inbounds A.data[row, col] = x
961✔
313
    end
1,232✔
314
    A
114✔
315
end
316

317
Base._reverse(A::UpperOrUnitUpperTriangular, dims::Integer) = reverse!(Matrix(A); dims)
8✔
318
Base._reverse(A::UpperTriangular, ::Colon) = LowerTriangular(reverse(A.data))
1✔
319
Base._reverse(A::UnitUpperTriangular, ::Colon) = UnitLowerTriangular(reverse(A.data))
1✔
320
Base._reverse(A::LowerOrUnitLowerTriangular, dims) = reverse!(Matrix(A); dims)
8✔
321
Base._reverse(A::LowerTriangular, ::Colon) = UpperTriangular(reverse(A.data))
1✔
322
Base._reverse(A::UnitLowerTriangular, ::Colon) = UnitUpperTriangular(reverse(A.data))
1✔
323

324
## structured matrix methods ##
325
function Base.replace_in_print_matrix(A::Union{UpperTriangular,UnitUpperTriangular},
158,562✔
326
                                      i::Integer, j::Integer, s::AbstractString)
327
    return i <= j ? s : Base.replace_with_centered_mark(s)
158,562✔
328
end
329
function Base.replace_in_print_matrix(A::Union{LowerTriangular,UnitLowerTriangular},
156,884✔
330
                                      i::Integer, j::Integer, s::AbstractString)
331
    return i >= j ? s : Base.replace_with_centered_mark(s)
156,884✔
332
end
333

334
istril(A::UnitLowerTriangular, k::Integer=0) = k >= 0
137✔
335
istriu(A::UnitUpperTriangular, k::Integer=0) = k <= 0
5,815✔
336
Base.@constprop :aggressive function istril(A::LowerTriangular, k::Integer=0)
116✔
337
    k >= 0 && return true
155✔
338
    return _istril(A, k)
45✔
339
end
340
@inline function _istril(A::LowerTriangular, k)
341
    P = parent(A)
9✔
342
    m = size(A, 1)
9✔
343
    for j in max(1, k + 2):m
9✔
344
        all(iszero, view(P, j:min(j - k - 1, m), j)) || return false
81✔
345
    end
153✔
346
    return true
9✔
347
end
348
Base.@constprop :aggressive function istriu(A::UpperTriangular, k::Integer=0)
123✔
349
    k <= 0 && return true
14,160✔
350
    return _istriu(A, k)
45✔
351
end
352
@inline function _istriu(A::UpperTriangular, k)
353
    P = parent(A)
9✔
354
    m = size(A, 1)
9✔
355
    for j in 1:min(m, m + k - 1)
9✔
356
        all(iszero, view(P, max(1, j - k + 1):j, j)) || return false
81✔
357
    end
153✔
358
    return true
9✔
359
end
360
istril(A::Adjoint, k::Integer=0) = istriu(A.parent, -k)
6,836✔
361
istril(A::Transpose, k::Integer=0) = istriu(A.parent, -k)
1,487✔
362
istriu(A::Adjoint, k::Integer=0) = istril(A.parent, -k)
6,847✔
363
istriu(A::Transpose, k::Integer=0) = istril(A.parent, -k)
1,487✔
364

365
function tril!(A::UpperTriangular{T}, k::Integer=0) where {T}
35✔
366
    n = size(A,1)
35✔
367
    if k < 0
35✔
368
        fill!(A.data, zero(T))
1,134✔
369
        return A
14✔
370
    elseif k == 0
21✔
371
        for j in 1:n, i in 1:j-1
7✔
372
            A.data[i,j] = zero(T)
252✔
373
        end
308✔
374
        return A
7✔
375
    else
376
        return UpperTriangular(tril!(A.data,k))
14✔
377
    end
378
end
379
function triu!(A::UpperTriangular, k::Integer=0)
132✔
380
    n = size(A,1)
194✔
381
    if k > 0
132✔
382
        for j in 1:n, i in max(1,j-k+1):j
28✔
383
            A.data[i,j] = zero(eltype(A))
756✔
384
        end
980✔
385
    end
386
    return A
132✔
387
end
388

389
function tril!(A::UnitUpperTriangular{T}, k::Integer=0) where {T}
35✔
390
    n = size(A,1)
35✔
391
    if k < 0
35✔
392
        fill!(A.data, zero(T))
1,134✔
393
        return UpperTriangular(A.data)
14✔
394
    elseif k == 0
21✔
395
        fill!(A.data, zero(T))
567✔
396
        for i in diagind(A)
14✔
397
            A.data[i] = oneunit(T)
63✔
398
        end
119✔
399
        return UpperTriangular(A.data)
7✔
400
    else
401
        for i in diagind(A)
28✔
402
            A.data[i] = oneunit(T)
126✔
403
        end
238✔
404
        return UpperTriangular(tril!(A.data,k))
14✔
405
    end
406
end
407

408
function triu!(A::UnitUpperTriangular, k::Integer=0)
10✔
409
    for i in diagind(A)
70✔
410
        A.data[i] = oneunit(eltype(A))
315✔
411
    end
595✔
412
    return triu!(UpperTriangular(A.data), k)
35✔
413
end
414

415
function triu!(A::LowerTriangular{T}, k::Integer=0) where {T}
35✔
416
    n = size(A,1)
35✔
417
    if k > 0
35✔
418
        fill!(A.data, zero(T))
1,134✔
419
        return A
14✔
420
    elseif k == 0
21✔
421
        for j in 1:n, i in j+1:n
14✔
422
            A.data[i,j] = zero(T)
252✔
423
        end
308✔
424
        return A
7✔
425
    else
426
        return LowerTriangular(triu!(A.data, k))
14✔
427
    end
428
end
429

430
function tril!(A::LowerTriangular, k::Integer=0)
70✔
431
    n = size(A,1)
70✔
432
    if k < 0
70✔
433
        for j in 1:n, i in j:min(j-k-1,n)
28✔
434
            A.data[i, j] = zero(eltype(A))
756✔
435
        end
980✔
436
    end
437
    A
70✔
438
end
439

440
function triu!(A::UnitLowerTriangular{T}, k::Integer=0) where T
35✔
441
    n = size(A,1)
35✔
442
    if k > 0
35✔
443
        fill!(A.data, zero(T))
1,134✔
444
        return LowerTriangular(A.data)
14✔
445
    elseif k == 0
21✔
446
        fill!(A.data, zero(T))
567✔
447
        for i in diagind(A)
14✔
448
            A.data[i] = oneunit(T)
63✔
449
        end
119✔
450
        return LowerTriangular(A.data)
7✔
451
    else
452
        for i in diagind(A)
28✔
453
            A.data[i] = oneunit(T)
126✔
454
        end
238✔
455
        return LowerTriangular(triu!(A.data, k))
14✔
456
    end
457
end
458

459
function tril!(A::UnitLowerTriangular, k::Integer=0)
10✔
460
    for i in diagind(A)
70✔
461
        A.data[i] = oneunit(eltype(A))
315✔
462
    end
595✔
463
    return tril!(LowerTriangular(A.data), k)
35✔
464
end
465

466
adjoint(A::LowerTriangular) = UpperTriangular(adjoint(A.data))
4,632✔
467
adjoint(A::UpperTriangular) = LowerTriangular(adjoint(A.data))
4,903✔
468
adjoint(A::UnitLowerTriangular) = UnitUpperTriangular(adjoint(A.data))
4,108✔
469
adjoint(A::UnitUpperTriangular) = UnitLowerTriangular(adjoint(A.data))
4,014✔
470
transpose(A::LowerTriangular) = UpperTriangular(transpose(A.data))
4,195✔
471
transpose(A::UpperTriangular) = LowerTriangular(transpose(A.data))
4,182✔
472
transpose(A::UnitLowerTriangular) = UnitUpperTriangular(transpose(A.data))
4,350✔
473
transpose(A::UnitUpperTriangular) = UnitLowerTriangular(transpose(A.data))
3,993✔
474

475
transpose!(A::LowerTriangular) = UpperTriangular(copytri!(A.data, 'L', false, true))
22✔
476
transpose!(A::UnitLowerTriangular) = UnitUpperTriangular(copytri!(A.data, 'L', false, false))
22✔
477
transpose!(A::UpperTriangular) = LowerTriangular(copytri!(A.data, 'U', false, true))
22✔
478
transpose!(A::UnitUpperTriangular) = UnitLowerTriangular(copytri!(A.data, 'U', false, false))
22✔
479
adjoint!(A::LowerTriangular) = UpperTriangular(copytri!(A.data, 'L' , true, true))
22✔
480
adjoint!(A::UnitLowerTriangular) = UnitUpperTriangular(copytri!(A.data, 'L' , true, false))
22✔
481
adjoint!(A::UpperTriangular) = LowerTriangular(copytri!(A.data, 'U' , true, true))
22✔
482
adjoint!(A::UnitUpperTriangular) = UnitLowerTriangular(copytri!(A.data, 'U' , true, false))
22✔
483

484
diag(A::UpperOrLowerTriangular) = diag(A.data)
275✔
485
diag(A::Union{UnitLowerTriangular, UnitUpperTriangular}) = fill(oneunit(eltype(A)), size(A,1))
484✔
486

487
# Unary operations
488
-(A::LowerTriangular) = LowerTriangular(-A.data)
1✔
489
-(A::UpperTriangular) = UpperTriangular(-A.data)
1✔
490
function -(A::UnitLowerTriangular)
11✔
491
    Adata = A.data
11✔
492
    Anew = similar(Adata) # must be mutable, even if Adata is not
22✔
493
    @. Anew = -Adata
22✔
494
    for i = 1:size(A, 1)
11✔
495
        Anew[i, i] = -A[i, i]
77✔
496
    end
143✔
497
    LowerTriangular(Anew)
11✔
498
end
499
function -(A::UnitUpperTriangular)
11✔
500
    Adata = A.data
11✔
501
    Anew = similar(Adata) # must be mutable, even if Adata is not
22✔
502
    @. Anew = -Adata
22✔
503
    for i = 1:size(A, 1)
11✔
504
        Anew[i, i] = -A[i, i]
77✔
505
    end
143✔
506
    UpperTriangular(Anew)
11✔
507
end
508

509
# use broadcasting if the parents are strided, where we loop only over the triangular part
510
for TM in (:LowerTriangular, :UpperTriangular)
511
    @eval -(A::$TM{<:Any, <:StridedMaybeAdjOrTransMat}) = broadcast(-, A)
205✔
512
end
513

514
tr(A::UpperOrLowerTriangular) = tr(A.data)
14✔
515
tr(A::Union{UnitLowerTriangular, UnitUpperTriangular}) = size(A, 1) * oneunit(eltype(A))
14✔
516

517
for T in (:UpperOrUnitUpperTriangular, :LowerOrUnitLowerTriangular)
518
    @eval @propagate_inbounds function copyto!(dest::$T, U::$T)
166✔
519
        if axes(dest) != axes(U)
12,584✔
520
            @invoke copyto!(dest::AbstractArray, U::AbstractArray)
20✔
521
        else
522
            _copyto!(dest, U)
6,287✔
523
        end
524
        return dest
6,285✔
525
    end
526
end
527

528
# copy and scale
529
for (T, UT) in ((:UpperTriangular, :UnitUpperTriangular), (:LowerTriangular, :UnitLowerTriangular))
530
    @eval @inline function _copyto!(A::$T, B::$T)
531
        @boundscheck checkbounds(A, axes(B)...)
5,717✔
532
        copytrito!(parent(A), parent(B), uplo_char(A))
5,717✔
533
        return A
5,717✔
534
    end
535
    @eval @inline function _copyto!(A::$UT, B::$T)
536
        for dind in diagind(A, IndexStyle(A))
8✔
537
            if A[dind] != B[dind]
16✔
538
                throw_nononeerror(typeof(A), B[dind], Tuple(dind)...)
2✔
539
            end
540
        end
26✔
541
        _copyto!($T(parent(A)), B)
2✔
542
        return A
2✔
543
    end
544
end
545
@inline function _copyto!(A::UpperOrUnitUpperTriangular, B::UnitUpperTriangular)
546
    @boundscheck checkbounds(A, axes(B)...)
296✔
547
    n = size(B,1)
296✔
548
    B2 = Base.unalias(A, B)
296✔
549
    for j = 1:n
296✔
550
        for i = 1:j-1
2,518✔
551
            @inbounds parent(A)[i,j] = parent(B2)[i,j]
9,699✔
552
        end
17,176✔
553
        if A isa UpperTriangular # copy diagonal
2,518✔
554
            @inbounds parent(A)[j,j] = B2[j,j]
12✔
555
        end
556
    end
4,740✔
557
    return A
296✔
558
end
559
@inline function _copyto!(A::LowerOrUnitLowerTriangular, B::UnitLowerTriangular)
560
    @boundscheck checkbounds(A, axes(B)...)
272✔
561
    n = size(B,1)
272✔
562
    B2 = Base.unalias(A, B)
272✔
563
    for j = 1:n
272✔
564
        if A isa LowerTriangular # copy diagonal
2,302✔
565
            @inbounds parent(A)[j,j] = B2[j,j]
12✔
566
        end
567
        for i = j+1:n
2,574✔
568
            @inbounds parent(A)[i,j] = parent(B2)[i,j]
8,835✔
569
        end
15,640✔
570
    end
4,332✔
571
    return A
272✔
572
end
573

574
_triangularize!(::UpperOrUnitUpperTriangular) = triu!
×
575
_triangularize!(::LowerOrUnitLowerTriangular) = tril!
×
576

577
@propagate_inbounds function copyto!(dest::StridedMatrix, U::UpperOrLowerTriangular)
8✔
578
    if axes(dest) != axes(U)
30,772✔
579
        @invoke copyto!(dest::StridedMatrix, U::AbstractArray)
×
580
    else
581
        _copyto!(dest, U)
15,386✔
582
    end
583
    return dest
15,386✔
584
end
585
@propagate_inbounds function _copyto!(dest::StridedMatrix, U::UpperOrLowerTriangular)
586
    copytrito!(dest, parent(U), U isa UpperOrUnitUpperTriangular ? 'U' : 'L')
1,259✔
587
    copytrito!(dest, U, U isa UpperOrUnitUpperTriangular ? 'L' : 'U')
1,259✔
588
    return dest
1,259✔
589
end
590
@propagate_inbounds function _copyto!(dest::StridedMatrix, U::UpperOrLowerTriangular{<:Any, <:StridedMatrix})
591
    U2 = Base.unalias(dest, U)
14,127✔
592
    copyto_unaliased!(dest, U2)
89,904✔
593
    return dest
14,127✔
594
end
595
# for strided matrices, we explicitly loop over the arrays to improve cache locality
596
# This fuses the copytrito! for the two halves
597
@inline function copyto_unaliased!(dest::StridedMatrix, U::UpperOrUnitUpperTriangular{<:Any, <:StridedMatrix})
598
    @boundscheck checkbounds(dest, axes(U)...)
10,405✔
599
    isunit = U isa UnitUpperTriangular
10,405✔
600
    for col in axes(dest,2)
10,405✔
601
        for row in 1:col-isunit
60,644✔
602
            @inbounds dest[row,col] = U.data[row,col]
257,740✔
603
        end
456,288✔
604
        for row in col+!isunit:size(U,1)
69,591✔
605
            @inbounds dest[row,col] = U[row,col]
216,562✔
606
        end
381,427✔
607
    end
110,889✔
608
    return dest
10,405✔
609
end
610
@inline function copyto_unaliased!(dest::StridedMatrix, L::LowerOrUnitLowerTriangular{<:Any, <:StridedMatrix})
611
    @boundscheck checkbounds(dest, axes(L)...)
3,722✔
612
    isunit = L isa UnitLowerTriangular
3,722✔
613
    for col in axes(dest,2)
3,722✔
614
        for row in 1:col-!isunit
29,248✔
615
            @inbounds dest[row,col] = L[row,col]
126,910✔
616
        end
227,260✔
617
        for row in col+isunit:size(L,1)
30,276✔
618
            @inbounds dest[row,col] = L.data[row,col]
139,956✔
619
        end
251,692✔
620
    end
54,780✔
621
    return dest
3,722✔
622
end
623

624
@inline _rscale_add!(A::AbstractTriangular, B::AbstractTriangular, C::Number, alpha::Number, beta::Number) =
620✔
625
    @stable_muladdmul _triscale!(A, B, C, MulAddMul(alpha, beta))
626
@inline _lscale_add!(A::AbstractTriangular, B::Number, C::AbstractTriangular, alpha::Number, beta::Number) =
26✔
627
    @stable_muladdmul _triscale!(A, B, C, MulAddMul(alpha, beta))
628

629
function checksize1(A, B)
630
    szA, szB = size(A), size(B)
2,118✔
631
    szA == szB || throw(DimensionMismatch(lazy"size of A, $szA, does not match size of B, $szB"))
2,126✔
632
    checksquare(B)
2,110✔
633
end
634

635
function _triscale!(A::UpperTriangular, B::UpperTriangular, c::Number, _add)
595✔
636
    n = checksize1(A, B)
835✔
637
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
834✔
638
    for j = 1:n
834✔
639
        for i = 1:j
3,200✔
640
            @inbounds _modify!(_add, B.data[i,j] * c, A.data, (i,j))
9,423✔
641
        end
15,646✔
642
    end
5,566✔
643
    return A
834✔
644
end
645
function _triscale!(A::UpperTriangular, c::Number, B::UpperTriangular, _add)
624✔
646
    n = checksize1(A, B)
663✔
647
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
662✔
648
    for j = 1:n
662✔
649
        for i = 1:j
2,879✔
650
            @inbounds _modify!(_add, c * B.data[i,j], A.data, (i,j))
10,432✔
651
        end
17,985✔
652
    end
5,096✔
653
    return A
662✔
654
end
655
function _triscale!(A::UpperOrUnitUpperTriangular, B::UnitUpperTriangular, c::Number, _add)
12✔
656
    n = checksize1(A, B)
14✔
657
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
13✔
658
    for j = 1:n
9✔
659
        @inbounds _modify!(_add, c, A, (j,j))
69✔
660
        for i = 1:(j - 1)
69✔
661
            @inbounds _modify!(_add, B.data[i,j] * c, A.data, (i,j))
258✔
662
        end
456✔
663
    end
129✔
664
    return A
9✔
665
end
666
function _triscale!(A::UpperOrUnitUpperTriangular, c::Number, B::UnitUpperTriangular, _add)
12✔
667
    n = checksize1(A, B)
12✔
668
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
11✔
669
    for j = 1:n
7✔
670
        @inbounds _modify!(_add, c, A, (j,j))
63✔
671
        for i = 1:(j - 1)
63✔
672
            @inbounds _modify!(_add, c * B.data[i,j], A.data, (i,j))
252✔
673
        end
448✔
674
    end
119✔
675
    return A
7✔
676
end
677
function _triscale!(A::LowerTriangular, B::LowerTriangular, c::Number, _add)
1✔
678
    n = checksize1(A, B)
96✔
679
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
95✔
680
    for j = 1:n
95✔
681
        for i = j:n
615✔
682
            @inbounds _modify!(_add, B.data[i,j] * c, A.data, (i,j))
2,795✔
683
        end
4,975✔
684
    end
1,135✔
685
    return A
95✔
686
end
687
function _triscale!(A::LowerTriangular, c::Number, B::LowerTriangular, _add)
277✔
688
    n = checksize1(A, B)
316✔
689
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
315✔
690
    for j = 1:n
315✔
691
        for i = j:n
1,611✔
692
            @inbounds _modify!(_add, c * B.data[i,j], A.data, (i,j))
6,323✔
693
        end
11,035✔
694
    end
2,907✔
695
    return A
315✔
696
end
697
function _triscale!(A::LowerOrUnitLowerTriangular, B::UnitLowerTriangular, c::Number, _add)
12✔
698
    n = checksize1(A, B)
14✔
699
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
13✔
700
    for j = 1:n
9✔
701
        @inbounds _modify!(_add, c, A, (j,j))
69✔
702
        for i = (j + 1):n
78✔
703
            @inbounds _modify!(_add, B.data[i,j] * c, A.data, (i,j))
258✔
704
        end
456✔
705
    end
129✔
706
    return A
9✔
707
end
708
function _triscale!(A::LowerOrUnitLowerTriangular, c::Number, B::UnitLowerTriangular, _add)
12✔
709
    n = checksize1(A, B)
12✔
710
    iszero(_add.alpha) && return _rmul_or_fill!(A, _add.beta)
11✔
711
    for j = 1:n
7✔
712
        @inbounds _modify!(_add, c, A, (j,j))
63✔
713
        for i = (j + 1):n
70✔
714
            @inbounds _modify!(_add, c * B.data[i,j], A.data, (i,j))
252✔
715
        end
448✔
716
    end
119✔
717
    return A
7✔
718
end
719

720
function _trirdiv!(A::UpperTriangular, B::UpperOrUnitUpperTriangular, c::Number)
94✔
721
    n = checksize1(A, B)
94✔
722
    for j in 1:n
94✔
723
        for i in 1:j
405✔
724
            @inbounds A[i, j] = B[i, j] / c
1,666✔
725
        end
2,927✔
726
    end
716✔
727
    return A
94✔
728
end
729
function _trirdiv!(A::LowerTriangular, B::LowerOrUnitLowerTriangular, c::Number)
22✔
730
    n = checksize1(A, B)
22✔
731
    for j in 1:n
22✔
732
        for i in j:n
142✔
733
            @inbounds A[i, j] = B[i, j] / c
654✔
734
        end
1,166✔
735
    end
262✔
736
    return A
22✔
737
end
738
function _trildiv!(A::UpperTriangular, c::Number, B::UpperOrUnitUpperTriangular)
20✔
739
    n = checksize1(A, B)
20✔
740
    for j in 1:n
20✔
741
        for i in 1:j
138✔
742
            @inbounds A[i, j] = c \ B[i, j]
648✔
743
        end
1,158✔
744
    end
256✔
745
    return A
20✔
746
end
747
function _trildiv!(A::LowerTriangular, c::Number, B::LowerOrUnitLowerTriangular)
20✔
748
    n = checksize1(A, B)
20✔
749
    for j in 1:n
20✔
750
        for i in j:n
138✔
751
            @inbounds A[i, j] = c \ B[i, j]
648✔
752
        end
1,158✔
753
    end
256✔
754
    return A
20✔
755
end
756

757
rmul!(A::UpperOrLowerTriangular, c::Number) = @inline _triscale!(A, A, c, MulAddMul())
339✔
758
lmul!(c::Number, A::UpperOrLowerTriangular) = @inline _triscale!(A, c, A, MulAddMul())
78✔
759

760
function dot(x::AbstractVector, A::UpperTriangular, y::AbstractVector)
42✔
761
    require_one_based_indexing(x, y)
42✔
762
    m = size(A, 1)
42✔
763
    (length(x) == m == length(y)) || throw(DimensionMismatch())
42✔
764
    if iszero(m)
42✔
765
        return dot(zero(eltype(x)), zero(eltype(A)), zero(eltype(y)))
×
766
    end
767
    x₁ = x[1]
42✔
768
    r = dot(x₁, A[1,1], y[1])
42✔
769
    @inbounds for j in 2:m
42✔
770
        yj = y[j]
336✔
771
        if !iszero(yj)
336✔
772
            temp = adjoint(A[1,j]) * x₁
336✔
773
            @simd for i in 2:j
336✔
774
                temp += adjoint(A[i,j]) * x[i]
1,512✔
775
            end
776
            r += dot(temp, yj)
336✔
777
        end
778
    end
630✔
779
    return r
42✔
780
end
781
function dot(x::AbstractVector, A::UnitUpperTriangular, y::AbstractVector)
42✔
782
    require_one_based_indexing(x, y)
42✔
783
    m = size(A, 1)
42✔
784
    (length(x) == m == length(y)) || throw(DimensionMismatch())
42✔
785
    if iszero(m)
42✔
786
        return dot(zero(eltype(x)), zero(eltype(A)), zero(eltype(y)))
×
787
    end
788
    x₁ = first(x)
42✔
789
    r = dot(x₁, y[1])
42✔
790
    @inbounds for j in 2:m
42✔
791
        yj = y[j]
336✔
792
        if !iszero(yj)
336✔
793
            temp = adjoint(A[1,j]) * x₁
336✔
794
            @simd for i in 2:j-1
336✔
795
                temp += adjoint(A[i,j]) * x[i]
1,176✔
796
            end
797
            r += dot(temp, yj)
420✔
798
            r += dot(x[j], yj)
336✔
799
        end
800
    end
630✔
801
    return r
42✔
802
end
803
function dot(x::AbstractVector, A::LowerTriangular, y::AbstractVector)
42✔
804
    require_one_based_indexing(x, y)
42✔
805
    m = size(A, 1)
42✔
806
    (length(x) == m == length(y)) || throw(DimensionMismatch())
42✔
807
    if iszero(m)
42✔
808
        return dot(zero(eltype(x)), zero(eltype(A)), zero(eltype(y)))
×
809
    end
810
    r = zero(typeof(dot(first(x), first(A), first(y))))
42✔
811
    @inbounds for j in 1:m
42✔
812
        yj = y[j]
378✔
813
        if !iszero(yj)
378✔
814
            temp = adjoint(A[j,j]) * x[j]
378✔
815
            @simd for i in j+1:m
420✔
816
                temp += adjoint(A[i,j]) * x[i]
1,512✔
817
            end
818
            r += dot(temp, yj)
378✔
819
        end
820
    end
714✔
821
    return r
42✔
822
end
823
function dot(x::AbstractVector, A::UnitLowerTriangular, y::AbstractVector)
42✔
824
    require_one_based_indexing(x, y)
42✔
825
    m = size(A, 1)
42✔
826
    (length(x) == m == length(y)) || throw(DimensionMismatch())
42✔
827
    if iszero(m)
42✔
828
        return dot(zero(eltype(x)), zero(eltype(A)), zero(eltype(y)))
×
829
    end
830
    r = zero(typeof(dot(first(x), first(y))))
42✔
831
    @inbounds for j in 1:m
42✔
832
        yj = y[j]
378✔
833
        if !iszero(yj)
378✔
834
            temp = x[j]
378✔
835
            @simd for i in j+1:m
420✔
836
                temp += adjoint(A[i,j]) * x[i]
1,876✔
837
            end
838
            r += dot(temp, yj)
491✔
839
        end
840
    end
714✔
841
    return r
42✔
842
end
843

844
fillstored!(A::LowerTriangular, x)     = (fillband!(A.data, x, 1-size(A,1), 0); A)
1✔
845
fillstored!(A::UnitLowerTriangular, x) = (fillband!(A.data, x, 1-size(A,1), -1); A)
1✔
846
fillstored!(A::UpperTriangular, x)     = (fillband!(A.data, x, 0, size(A,2)-1); A)
1✔
847
fillstored!(A::UnitUpperTriangular, x) = (fillband!(A.data, x, 1, size(A,2)-1); A)
×
848

849
# Binary operations
850
# use broadcasting if the parents are strided, where we loop only over the triangular part
851
function +(A::UpperTriangular, B::UpperTriangular)
245✔
852
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
271✔
853
    UpperTriangular(A.data + B.data)
11✔
854
end
855
function +(A::LowerTriangular, B::LowerTriangular)
66✔
856
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
72✔
857
    LowerTriangular(A.data + B.data)
7✔
858
end
859
function +(A::UpperTriangular, B::UnitUpperTriangular)
55✔
860
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
55✔
861
    UpperTriangular(A.data + triu(B.data, 1) + I)
1✔
862
end
863
function +(A::LowerTriangular, B::UnitLowerTriangular)
51✔
864
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
51✔
865
    LowerTriangular(A.data + tril(B.data, -1) + I)
1✔
866
end
867
function +(A::UnitUpperTriangular, B::UpperTriangular)
55✔
868
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
55✔
869
    UpperTriangular(triu(A.data, 1) + B.data + I)
1✔
870
end
871
function +(A::UnitLowerTriangular, B::LowerTriangular)
50✔
872
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
50✔
873
    LowerTriangular(tril(A.data, -1) + B.data + I)
1✔
874
end
875
function +(A::UnitUpperTriangular, B::UnitUpperTriangular)
58✔
876
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
58✔
877
    UpperTriangular(triu(A.data, 1) + triu(B.data, 1) + 2I)
1✔
878
end
879
function +(A::UnitLowerTriangular, B::UnitLowerTriangular)
50✔
880
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .+ B
50✔
881
    LowerTriangular(tril(A.data, -1) + tril(B.data, -1) + 2I)
1✔
882
end
883
+(A::AbstractTriangular, B::AbstractTriangular) = copyto!(similar(parent(A)), A) + copyto!(similar(parent(B)), B)
406✔
884

885
function -(A::UpperTriangular, B::UpperTriangular)
338✔
886
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
380✔
887
    UpperTriangular(A.data - B.data)
13✔
888
end
889
function -(A::LowerTriangular, B::LowerTriangular)
267✔
890
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
274✔
891
    LowerTriangular(A.data - B.data)
3✔
892
end
893
function -(A::UpperTriangular, B::UnitUpperTriangular)
51✔
894
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
51✔
895
    UpperTriangular(A.data - triu(B.data, 1) - I)
1✔
896
end
897
function -(A::LowerTriangular, B::UnitLowerTriangular)
50✔
898
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
50✔
899
    LowerTriangular(A.data - tril(B.data, -1) - I)
1✔
900
end
901
function -(A::UnitUpperTriangular, B::UpperTriangular)
51✔
902
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
51✔
903
    UpperTriangular(triu(A.data, 1) - B.data + I)
1✔
904
end
905
function -(A::UnitLowerTriangular, B::LowerTriangular)
50✔
906
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
50✔
907
    LowerTriangular(tril(A.data, -1) - B.data + I)
1✔
908
end
909
function -(A::UnitUpperTriangular, B::UnitUpperTriangular)
52✔
910
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
58✔
911
    UpperTriangular(triu(A.data, 1) - triu(B.data, 1))
1✔
912
end
913
function -(A::UnitLowerTriangular, B::UnitLowerTriangular)
52✔
914
    (parent(A) isa StridedMatrix || parent(B) isa StridedMatrix) && return A .- B
57✔
915
    LowerTriangular(tril(A.data, -1) - tril(B.data, -1))
1✔
916
end
917
-(A::AbstractTriangular, B::AbstractTriangular) = copyto!(similar(parent(A)), A) - copyto!(similar(parent(B)), B)
396✔
918

919
function kron(A::UpperTriangular{T,<:StridedMaybeAdjOrTransMat}, B::UpperTriangular{S,<:StridedMaybeAdjOrTransMat}) where {T,S}
57✔
920
    C = UpperTriangular(Matrix{promote_op(*, T, S)}(undef, _kronsize(A, B)))
57✔
921
    return kron!(C, A, B)
57✔
922
end
923
function kron(A::LowerTriangular{T,<:StridedMaybeAdjOrTransMat}, B::LowerTriangular{S,<:StridedMaybeAdjOrTransMat}) where {T,S}
57✔
924
    C = LowerTriangular(Matrix{promote_op(*, T, S)}(undef, _kronsize(A, B)))
57✔
925
    return kron!(C, A, B)
57✔
926
end
927

928
function kron!(C::UpperTriangular{<:Any,<:StridedMaybeAdjOrTransMat}, A::UpperTriangular{<:Any,<:StridedMaybeAdjOrTransMat}, B::UpperTriangular{<:Any,<:StridedMaybeAdjOrTransMat})
929
    size(C) == _kronsize(A, B) || throw(DimensionMismatch("kron!"))
57✔
930
    _triukron!(C.data, A.data, B.data)
57✔
931
    return C
57✔
932
end
933
function kron!(C::LowerTriangular{<:Any,<:StridedMaybeAdjOrTransMat}, A::LowerTriangular{<:Any,<:StridedMaybeAdjOrTransMat}, B::LowerTriangular{<:Any,<:StridedMaybeAdjOrTransMat})
934
    size(C) == _kronsize(A, B) || throw(DimensionMismatch("kron!"))
57✔
935
    _trilkron!(C.data, A.data, B.data)
57✔
936
    return C
57✔
937
end
938

939
function _triukron!(C, A, B)
57✔
940
    n_A = size(A, 1)
57✔
941
    n_B = size(B, 1)
57✔
942
    @inbounds for j = 1:n_A
57✔
943
        jnB = (j - 1) * n_B
457✔
944
        for i = 1:(j-1)
457✔
945
            Aij = A[i, j]
1,772✔
946
            inB = (i - 1) * n_B
1,772✔
947
            for l = 1:n_B
1,772✔
948
                for k = 1:l
15,892✔
949
                    C[inB+k, jnB+l] = Aij * B[k, l]
79,404✔
950
                end
142,916✔
951
                for k = 1:(l-1)
15,892✔
952
                    C[inB+l, jnB+k] = zero(C[inB+k, jnB+l])
63,512✔
953
                end
112,904✔
954
            end
30,012✔
955
        end
3,144✔
956
        Ajj = A[j, j]
457✔
957
        for l = 1:n_B
457✔
958
            for k = 1:l
4,001✔
959
                C[jnB+k, jnB+l] = Ajj * B[k, l]
19,893✔
960
            end
35,785✔
961
        end
7,545✔
962
    end
457✔
963
end
964

965
function _trilkron!(C, A, B)
57✔
966
    n_A = size(A, 1)
57✔
967
    n_B = size(B, 1)
57✔
968
    @inbounds for j = 1:n_A
57✔
969
        jnB = (j - 1) * n_B
457✔
970
        Ajj = A[j, j]
457✔
971
        for l = 1:n_B
457✔
972
            for k = l:n_B
4,001✔
973
                C[jnB+k, jnB+l] = Ajj * B[k, l]
19,893✔
974
            end
35,785✔
975
        end
7,545✔
976
        for i = (j+1):n_A
514✔
977
            Aij = A[i, j]
1,772✔
978
            inB = (i - 1) * n_B
1,772✔
979
            for l = 1:n_B
1,772✔
980
                for k = l:n_B
15,892✔
981
                    C[inB+k, jnB+l] = Aij * B[k, l]
79,404✔
982
                end
142,916✔
983
                for k = (l+1):n_B
17,664✔
984
                    C[inB+l, jnB+k] = zero(C[inB+k, jnB+l])
63,512✔
985
                end
112,904✔
986
            end
30,012✔
987
        end
3,144✔
988
    end
457✔
989
end
990

991
######################
992
# BlasFloat routines #
993
######################
994

995
# which triangle to use of the underlying data
996
uplo_char(::UpperOrUnitUpperTriangular) = 'U'
×
997
uplo_char(::LowerOrUnitLowerTriangular) = 'L'
×
998
uplo_char(::UpperOrUnitUpperTriangular{<:Any,<:AdjOrTrans}) = 'L'
×
999
uplo_char(::LowerOrUnitLowerTriangular{<:Any,<:AdjOrTrans}) = 'U'
×
1000
uplo_char(::UpperOrUnitUpperTriangular{<:Any,<:Adjoint{<:Any,<:Transpose}}) = 'U'
×
1001
uplo_char(::LowerOrUnitLowerTriangular{<:Any,<:Adjoint{<:Any,<:Transpose}}) = 'L'
×
1002
uplo_char(::UpperOrUnitUpperTriangular{<:Any,<:Transpose{<:Any,<:Adjoint}}) = 'U'
×
1003
uplo_char(::LowerOrUnitLowerTriangular{<:Any,<:Transpose{<:Any,<:Adjoint}}) = 'L'
×
1004

1005
isunit_char(::UpperTriangular) = 'N'
×
1006
isunit_char(::UnitUpperTriangular) = 'U'
×
1007
isunit_char(::LowerTriangular) = 'N'
×
1008
isunit_char(::UnitLowerTriangular) = 'U'
×
1009

1010
# generic fallback for AbstractTriangular matrices outside of the four subtypes provided here
1011
_trimul!(C::AbstractVecOrMat, A::AbstractTriangular, B::AbstractVector) =
×
1012
    lmul!(A, copyto!(C, B))
1013
_trimul!(C::AbstractMatrix, A::AbstractTriangular, B::AbstractMatrix) =
10✔
1014
    lmul!(A, copyto!(C, B))
1015
_trimul!(C::AbstractMatrix, A::AbstractMatrix, B::AbstractTriangular) =
4✔
1016
    rmul!(copyto!(C, A), B)
1017
_trimul!(C::AbstractMatrix, A::AbstractTriangular, B::AbstractTriangular) =
×
1018
    lmul!(A, copyto!(C, B))
1019
# redirect for UpperOrLowerTriangular
1020
_trimul!(C::AbstractVecOrMat, A::UpperOrLowerTriangular, B::AbstractVector) =
3,481✔
1021
    generic_trimatmul!(C, uplo_char(A), isunit_char(A), wrapperop(parent(A)), _unwrap_at(parent(A)), B)
1022
_trimul!(C::AbstractMatrix, A::UpperOrLowerTriangular, B::AbstractMatrix) =
5,600✔
1023
    generic_trimatmul!(C, uplo_char(A), isunit_char(A), wrapperop(parent(A)), _unwrap_at(parent(A)), B)
1024
_trimul!(C::AbstractMatrix, A::AbstractMatrix, B::UpperOrLowerTriangular) =
6,433✔
1025
    generic_mattrimul!(C, uplo_char(B), isunit_char(B), wrapperop(parent(B)), A, _unwrap_at(parent(B)))
1026
_trimul!(C::AbstractMatrix, A::UpperOrLowerTriangular, B::UpperOrLowerTriangular) =
11,963✔
1027
    generic_trimatmul!(C, uplo_char(A), isunit_char(A), wrapperop(parent(A)), _unwrap_at(parent(A)), B)
1028
# disambiguation with AbstractTriangular
1029
_trimul!(C::AbstractMatrix, A::UpperOrLowerTriangular, B::AbstractTriangular) =
×
1030
    generic_trimatmul!(C, uplo_char(A), isunit_char(A), wrapperop(parent(A)), _unwrap_at(parent(A)), B)
1031
_trimul!(C::AbstractMatrix, A::AbstractTriangular, B::UpperOrLowerTriangular) =
×
1032
    generic_mattrimul!(C, uplo_char(B), isunit_char(B), wrapperop(parent(B)), A, _unwrap_at(parent(B)))
1033

1034
function lmul!(A::AbstractTriangular, B::AbstractVecOrMat)
556✔
1035
    if istriu(A)
804✔
1036
        _trimul!(B, UpperTriangular(A), B)
282✔
1037
    else
1038
        _trimul!(B, LowerTriangular(A), B)
522✔
1039
    end
1040
end
1041
function rmul!(A::AbstractMatrix, B::AbstractTriangular)
516✔
1042
    if istriu(B)
560✔
1043
        _trimul!(A, A, UpperTriangular(B))
280✔
1044
    else
1045
        _trimul!(A, A, LowerTriangular(B))
280✔
1046
    end
1047
end
1048

1049
for TC in (:AbstractVector, :AbstractMatrix)
1050
    @eval @inline function _mul!(C::$TC, A::AbstractTriangular, B::AbstractVector, alpha::Number, beta::Number)
1051
        if isone(alpha) && iszero(beta)
2,787✔
1052
            return _trimul!(C, A, B)
2,595✔
1053
        else
1054
            return @stable_muladdmul generic_matvecmul!(C, 'N', A, B, MulAddMul(alpha, beta))
268✔
1055
        end
1056
    end
1057
end
1058
for (TA, TB) in ((:AbstractTriangular, :AbstractMatrix),
1059
                    (:AbstractMatrix, :AbstractTriangular),
1060
                    (:AbstractTriangular, :AbstractTriangular)
1061
                )
1062
    @eval @inline function _mul!(C::AbstractMatrix, A::$TA, B::$TB, alpha::Number, beta::Number)
1063
        if isone(alpha) && iszero(beta)
21,999✔
1064
            return _trimul!(C, A, B)
21,868✔
1065
        else
1066
            return generic_matmatmul!(C, 'N', 'N', A, B, alpha, beta)
131✔
1067
        end
1068
    end
1069
end
1070

1071
ldiv!(C::AbstractVecOrMat, A::AbstractTriangular, B::AbstractVecOrMat) = _ldiv!(C, A, B)
15,786✔
1072
# generic fallback for AbstractTriangular, directs to 2-arg [l/r]div!
1073
_ldiv!(C::AbstractVecOrMat, A::AbstractTriangular, B::AbstractVecOrMat) =
×
1074
    ldiv!(A, copyto!(C, B))
1075
_rdiv!(C::AbstractMatrix, A::AbstractMatrix, B::AbstractTriangular) =
×
1076
    rdiv!(copyto!(C, A), B)
1077
# redirect for UpperOrLowerTriangular to generic_*div!
1078
_ldiv!(C::AbstractVecOrMat, A::UpperOrLowerTriangular, B::AbstractVecOrMat) =
28,394✔
1079
    generic_trimatdiv!(C, uplo_char(A), isunit_char(A), wrapperop(parent(A)), _unwrap_at(parent(A)), B)
1080
_rdiv!(C::AbstractMatrix, A::AbstractMatrix, B::UpperOrLowerTriangular) =
7,360✔
1081
    generic_mattridiv!(C, uplo_char(B), isunit_char(B), wrapperop(parent(B)), A, _unwrap_at(parent(B)))
1082

1083
function ldiv!(A::AbstractTriangular, B::AbstractVecOrMat)
6,026✔
1084
    if istriu(A)
12,742✔
1085
        _ldiv!(B, UpperTriangular(A), B)
7,113✔
1086
    else
1087
        _ldiv!(B, LowerTriangular(A), B)
5,495✔
1088
    end
1089
end
1090
function rdiv!(A::AbstractMatrix, B::AbstractTriangular)
104✔
1091
    if istriu(B)
1,956✔
1092
        _rdiv!(A, A, UpperTriangular(B))
1,844✔
1093
    else
1094
        _rdiv!(A, A, LowerTriangular(B))
112✔
1095
    end
1096
end
1097

1098
# preserve triangular structure in in-place multiplication/division
1099
for (cty, aty, bty) in ((:UpperTriangular, :UpperTriangular, :UpperTriangular),
1100
                        (:UpperTriangular, :UpperTriangular, :UnitUpperTriangular),
1101
                        (:UpperTriangular, :UnitUpperTriangular, :UpperTriangular),
1102
                        (:UnitUpperTriangular, :UnitUpperTriangular, :UnitUpperTriangular),
1103
                        (:LowerTriangular, :LowerTriangular, :LowerTriangular),
1104
                        (:LowerTriangular, :LowerTriangular, :UnitLowerTriangular),
1105
                        (:LowerTriangular, :UnitLowerTriangular, :LowerTriangular),
1106
                        (:UnitLowerTriangular, :UnitLowerTriangular, :UnitLowerTriangular))
1107
    @eval begin
1108
        function _trimul!(C::$cty, A::$aty, B::$bty)
1109
            _trimul!(parent(C), A, B)
1,302✔
1110
            return C
1,302✔
1111
        end
1112
        function _ldiv!(C::$cty, A::$aty, B::$bty)
1113
            _ldiv!(parent(C), A, B)
666✔
1114
            return C
666✔
1115
        end
1116
        function _rdiv!(C::$cty, A::$aty, B::$bty)
1117
            _rdiv!(parent(C), A, B)
1,705✔
1118
            return C
1,705✔
1119
        end
1120
    end
1121
end
1122

1123
for (t, uploc, isunitc) in ((:LowerTriangular, 'L', 'N'),
1124
                            (:UnitLowerTriangular, 'L', 'U'),
1125
                            (:UpperTriangular, 'U', 'N'),
1126
                            (:UnitUpperTriangular, 'U', 'U'))
1127
    @eval begin
1128
        # Matrix inverse
1129
        inv!(A::$t{T,S}) where {T<:BlasFloat,S<:StridedMatrix} =
16✔
1130
            $t{T,S}(LAPACK.trtri!($uploc, $isunitc, A.data))
1131

1132
        function inv(A::$t{T}) where {T}
173✔
1133
            S = typeof(inv(oneunit(T)))
533✔
1134
            if S <: BlasFloat || S === T # i.e. A is unitless
533✔
1135
                $t(ldiv!(convert(AbstractArray{S}, A), Matrix{S}(I, size(A))))
464✔
1136
            else
1137
                J = (one(T)*I)(size(A, 1))
650✔
1138
                $t(ldiv!(similar(A, S, size(A)), A, J))
69✔
1139
            end
1140
        end
1141

1142
        # Error bounds for triangular solve
1143
        errorbounds(A::$t{T,<:StridedMatrix}, X::StridedVecOrMat{T}, B::StridedVecOrMat{T}) where {T<:BlasFloat} =
160✔
1144
            LAPACK.trrfs!($uploc, 'N', $isunitc, A.data, B, X)
1145

1146
        # Condition numbers
1147
        function cond(A::$t{<:BlasFloat,<:StridedMatrix}, p::Real=2)
144✔
1148
            checksquare(A)
144✔
1149
            if p == 1
144✔
1150
                return inv(LAPACK.trcon!('O', $uploc, $isunitc, A.data))
64✔
1151
            elseif p == Inf
80✔
1152
                return inv(LAPACK.trcon!('I', $uploc, $isunitc, A.data))
64✔
1153
            else # use fallback
1154
                return cond(copyto!(similar(parent(A)), A), p)
16✔
1155
            end
1156
        end
1157
    end
1158
end
1159

1160
# multiplication
1161
generic_trimatmul!(c::StridedVector{T}, uploc, isunitc, tfun::Function, A::StridedMatrix{T}, b::AbstractVector{T}) where {T<:BlasFloat} =
952✔
1162
    BLAS.trmv!(uploc, tfun === identity ? 'N' : tfun === transpose ? 'T' : 'C', isunitc, A, c === b ? c : copyto!(c, b))
1163
generic_trimatmul!(C::StridedMatrix{T}, uploc, isunitc, tfun::Function, A::StridedMatrix{T}, B::AbstractMatrix{T}) where {T<:BlasFloat} =
6,273✔
1164
    BLAS.trmm!('L', uploc, tfun === identity ? 'N' : tfun === transpose ? 'T' : 'C', isunitc, one(T), A, C === B ? C : copyto!(C, B))
1165
generic_mattrimul!(C::StridedMatrix{T}, uploc, isunitc, tfun::Function, A::AbstractMatrix{T}, B::StridedMatrix{T}) where {T<:BlasFloat} =
1,664✔
1166
    BLAS.trmm!('R', uploc, tfun === identity ? 'N' : tfun === transpose ? 'T' : 'C', isunitc, one(T), B, C === A ? C : copyto!(C, A))
1167
# division
1168
generic_trimatdiv!(C::StridedVecOrMat{T}, uploc, isunitc, tfun::Function, A::StridedMatrix{T}, B::AbstractVecOrMat{T}) where {T<:BlasFloat} =
4,758✔
1169
    LAPACK.trtrs!(uploc, tfun === identity ? 'N' : tfun === transpose ? 'T' : 'C', isunitc, A, C === B ? C : copyto!(C, B))
1170
generic_mattridiv!(C::StridedMatrix{T}, uploc, isunitc, tfun::Function, A::AbstractMatrix{T}, B::StridedMatrix{T}) where {T<:BlasFloat} =
2,500✔
1171
    BLAS.trsm!('R', uploc, tfun === identity ? 'N' : tfun === transpose ? 'T' : 'C', isunitc, one(T), B, C === A ? C : copyto!(C, A))
1172

1173
errorbounds(A::AbstractTriangular{T}, X::AbstractVecOrMat{T}, B::AbstractVecOrMat{T}) where {T<:Union{BigFloat,Complex{BigFloat}}} =
×
1174
    error("not implemented yet! Please submit a pull request.")
1175
function errorbounds(A::AbstractTriangular{TA}, X::AbstractVecOrMat{TX}, B::AbstractVecOrMat{TB}) where {TA<:Number,TX<:Number,TB<:Number}
128✔
1176
    TAXB = promote_type(TA, TB, TX, Float32)
128✔
1177
    errorbounds(convert(AbstractMatrix{TAXB}, A), convert(AbstractArray{TAXB}, X), convert(AbstractArray{TAXB}, B))
128✔
1178
end
1179

1180
# Eigensystems
1181
## Notice that trecv works for quasi-triangular matrices and therefore the lower sub diagonal must be zeroed before calling the subroutine
1182
function eigvecs(A::UpperTriangular{<:BlasFloat,<:StridedMatrix})
1✔
1183
    LAPACK.trevc!('R', 'A', BlasInt[], triu!(A.data))
10✔
1184
end
1185
function eigvecs(A::UnitUpperTriangular{<:BlasFloat,<:StridedMatrix})
5✔
1186
    for i = 1:size(A, 1)
5✔
1187
        A.data[i,i] = 1
45✔
1188
    end
85✔
1189
    LAPACK.trevc!('R', 'A', BlasInt[], triu!(A.data))
5✔
1190
end
1191
function eigvecs(A::LowerTriangular{<:BlasFloat,<:StridedMatrix})
9✔
1192
    LAPACK.trevc!('L', 'A', BlasInt[], copy(tril!(A.data)'))
9✔
1193
end
1194
function eigvecs(A::UnitLowerTriangular{<:BlasFloat,<:StridedMatrix})
5✔
1195
    for i = 1:size(A, 1)
5✔
1196
        A.data[i,i] = 1
45✔
1197
    end
85✔
1198
    LAPACK.trevc!('L', 'A', BlasInt[], copy(tril!(A.data)'))
5✔
1199
end
1200

1201
####################
1202
# Generic routines #
1203
####################
1204

1205
for (t, unitt) in ((UpperTriangular, UnitUpperTriangular),
1206
                   (LowerTriangular, UnitLowerTriangular))
1207
    tstrided = t{<:Any, <:StridedMaybeAdjOrTransMat}
1208
    @eval begin
1209
        (*)(A::$t, x::Number) = $t(A.data*x)
4✔
1210
        function (*)(A::$tstrided, x::Number)
20✔
1211
            eltype_dest = promote_op(*, eltype(A), typeof(x))
76✔
1212
            dest = $t(similar(parent(A), eltype_dest))
76✔
1213
            _triscale!(dest, x, A, MulAddMul())
76✔
1214
        end
1215

1216
        function (*)(A::$unitt, x::Number)
26✔
1217
            B = $t(A.data)*x
26✔
1218
            for i = 1:size(A, 1)
26✔
1219
                B.data[i,i] = x
192✔
1220
            end
358✔
1221
            return B
26✔
1222
        end
1223

1224
        (*)(x::Number, A::$t) = $t(x*A.data)
×
1225
        function (*)(x::Number, A::$tstrided)
575✔
1226
            eltype_dest = promote_op(*, typeof(x), eltype(A))
823✔
1227
            dest = $t(similar(parent(A), eltype_dest))
823✔
1228
            _triscale!(dest, x, A, MulAddMul())
823✔
1229
        end
1230

1231
        function (*)(x::Number, A::$unitt)
42✔
1232
            B = x*$t(A.data)
42✔
1233
            for i = 1:size(A, 1)
42✔
1234
                B.data[i,i] = x
336✔
1235
            end
630✔
1236
            return B
42✔
1237
        end
1238

1239
        (/)(A::$t, x::Number) = $t(A.data/x)
2✔
1240
        function (/)(A::$tstrided, x::Number)
20✔
1241
            eltype_dest = promote_op(/, eltype(A), typeof(x))
116✔
1242
            dest = $t(similar(parent(A), eltype_dest))
116✔
1243
            _trirdiv!(dest, A,  x)
116✔
1244
        end
1245

1246
        function (/)(A::$unitt, x::Number)
20✔
1247
            B = $t(A.data)/x
20✔
1248
            invx = inv(x)
20✔
1249
            for i = 1:size(A, 1)
20✔
1250
                B.data[i,i] = invx
138✔
1251
            end
256✔
1252
            return B
20✔
1253
        end
1254

1255
        (\)(x::Number, A::$t) = $t(x\A.data)
×
1256
        function (\)(x::Number, A::$tstrided)
20✔
1257
            eltype_dest = promote_op(\, typeof(x), eltype(A))
40✔
1258
            dest = $t(similar(parent(A), eltype_dest))
40✔
1259
            _trildiv!(dest, x, A)
40✔
1260
        end
1261

1262
        function (\)(x::Number, A::$unitt)
20✔
1263
            B = x\$t(A.data)
20✔
1264
            invx = inv(x)
20✔
1265
            for i = 1:size(A, 1)
20✔
1266
                B.data[i,i] = invx
138✔
1267
            end
256✔
1268
            return B
20✔
1269
        end
1270
    end
1271
end
1272

1273
## Generic triangular multiplication
1274
function generic_trimatmul!(C::AbstractVecOrMat, uploc, isunitc, tfun::Function, A::AbstractMatrix, B::AbstractVecOrMat)
13,597✔
1275
    require_one_based_indexing(C, A, B)
13,597✔
1276
    m, n = size(B, 1), size(B, 2)
13,597✔
1277
    N = size(A, 1)
13,597✔
1278
    if m != N
13,597✔
1279
        throw(DimensionMismatch(lazy"right hand side B needs first dimension of size $(size(A,1)), has size $m"))
2,496✔
1280
    end
1281
    mc, nc = size(C, 1), size(C, 2)
11,101✔
1282
    if mc != N || nc != n
22,202✔
1283
        throw(DimensionMismatch(lazy"output has dimensions ($mc,$nc), should have ($N,$n)"))
×
1284
    end
1285
    oA = oneunit(eltype(A))
11,101✔
1286
    unit = isunitc == 'U'
11,101✔
1287
    @inbounds if uploc == 'U'
11,101✔
1288
        if tfun === identity
5,736✔
1289
            for j in 1:n
2,743✔
1290
                for i in 1:m
17,941✔
1291
                    Cij = (unit ? oA : A[i,i]) * B[i,j]
159,628✔
1292
                    for k in i + 1:m
176,869✔
1293
                        Cij += A[i,k] * B[k,j]
638,098✔
1294
                    end
1,132,409✔
1295
                    C[i,j] = Cij
158,928✔
1296
                end
299,915✔
1297
            end
33,139✔
1298
        else # tfun in (transpose, adjoint)
1299
            for j in 1:n
2,993✔
1300
                for i in m:-1:1
46,909✔
1301
                    Cij = (unit ? oA : tfun(A[i,i])) * B[i,j]
211,616✔
1302
                    for k in 1:i - 1
210,266✔
1303
                        Cij += tfun(A[k,i]) * B[k,j]
842,474✔
1304
                    end
1,492,737✔
1305
                    C[i,j] = Cij
210,266✔
1306
                end
397,077✔
1307
            end
43,917✔
1308
        end
1309
    else # uploc == 'L'
1310
        if tfun === identity
5,365✔
1311
            for j in 1:n
2,294✔
1312
                for i in m:-1:1
33,315✔
1313
                    Cij = (unit ? oA : A[i,i]) * B[i,j]
150,181✔
1314
                    for k in 1:i - 1
148,081✔
1315
                        Cij += A[i,k] * B[k,j]
596,054✔
1316
                    end
1,052,285✔
1317
                    C[i,j] = Cij
148,081✔
1318
                end
279,504✔
1319
            end
16,658✔
1320
        else # tfun in (transpose, adjoint)
1321
            for j in 1:n
3,071✔
1322
                for i in 1:m
23,719✔
1323
                    Cij = (unit ? oA : tfun(A[i,i])) * B[i,j]
214,302✔
1324
                    for k in i + 1:m
233,971✔
1325
                        Cij += tfun(A[k,i]) * B[k,j]
841,113✔
1326
                    end
1,479,493✔
1327
                    C[i,j] = Cij
210,252✔
1328
                end
396,785✔
1329
            end
23,719✔
1330
        end
1331
    end
1332
    return C
11,101✔
1333
end
1334
# conjugate cases
1335
function generic_trimatmul!(C::AbstractVecOrMat, uploc, isunitc, ::Function, xA::AdjOrTrans, B::AbstractVecOrMat)
22✔
1336
    A = parent(xA)
22✔
1337
    require_one_based_indexing(C, A, B)
22✔
1338
    m, n = size(B, 1), size(B, 2)
22✔
1339
    N = size(A, 1)
22✔
1340
    if m != N
22✔
1341
        throw(DimensionMismatch(lazy"right hand side B needs first dimension of size $(size(A,1)), has size $m"))
×
1342
    end
1343
    mc, nc = size(C, 1), size(C, 2)
22✔
1344
    if mc != N || nc != n
44✔
1345
        throw(DimensionMismatch(lazy"output has dimensions ($mc,$nc), should have ($N,$n)"))
×
1346
    end
1347
    oA = oneunit(eltype(A))
22✔
1348
    unit = isunitc == 'U'
22✔
1349
    @inbounds if uploc == 'U'
22✔
1350
        for j in 1:n
11✔
1351
            for i in 1:m
11✔
1352
                Cij = (unit ? oA : conj(A[i,i])) * B[i,j]
4,040✔
1353
                for k in i + 1:m
4,051✔
1354
                    Cij += conj(A[i,k]) * B[k,j]
1,998,147✔
1355
                end
3,992,265✔
1356
                C[i,j] = Cij
4,040✔
1357
            end
8,069✔
1358
        end
11✔
1359
    else # uploc == 'L'
1360
        for j in 1:n
11✔
1361
            for i in m:-1:1
20✔
1362
                Cij = (unit ? oA : conj(A[i,i])) * B[i,j]
4,040✔
1363
                for k in 1:i - 1
4,040✔
1364
                    Cij += conj(A[i,k]) * B[k,j]
1,998,147✔
1365
                end
3,992,265✔
1366
                C[i,j] = Cij
4,040✔
1367
            end
8,069✔
1368
        end
11✔
1369
    end
1370
    return C
22✔
1371
end
1372

1373
function generic_mattrimul!(C::AbstractMatrix, uploc, isunitc, tfun::Function, A::AbstractMatrix, B::AbstractMatrix)
4,769✔
1374
    require_one_based_indexing(C, A, B)
4,769✔
1375
    m, n = size(A, 1), size(A, 2)
4,769✔
1376
    N = size(B, 1)
4,769✔
1377
    if n != N
4,769✔
1378
        throw(DimensionMismatch(lazy"right hand side B needs first dimension of size $n, has size $N"))
2,496✔
1379
    end
1380
    mc, nc = size(C, 1), size(C, 2)
2,273✔
1381
    if mc != m || nc != N
4,546✔
1382
        throw(DimensionMismatch(lazy"output has dimensions ($mc,$nc), should have ($m,$N)"))
×
1383
    end
1384
    oB = oneunit(eltype(B))
2,273✔
1385
    unit = isunitc == 'U'
2,273✔
1386
    @inbounds if uploc == 'U'
2,273✔
1387
        if tfun === identity
1,079✔
1388
            for i in 1:m
444✔
1389
                for j in n:-1:1
7,428✔
1390
                    Cij = A[i,j] * (unit ? oB : B[j,j])
51,454✔
1391
                    for k in 1:j - 1
33,344✔
1392
                        Cij += A[i,k] * B[k,j]
137,210✔
1393
                    end
244,790✔
1394
                    C[i,j] = Cij
33,344✔
1395
                end
62,974✔
1396
            end
6,984✔
1397
        else # tfun in (transpose, adjoint)
1398
            for i in 1:m
635✔
1399
                for j in 1:n
5,613✔
1400
                    Cij = A[i,j] * (unit ? oB : tfun(B[j,j]))
70,761✔
1401
                    for k in j + 1:n
56,232✔
1402
                        Cij += A[i,k] * tfun(B[j,k])
205,749✔
1403
                    end
366,492✔
1404
                    C[i,j] = Cij
50,619✔
1405
                end
95,625✔
1406
            end
10,591✔
1407
        end
1408
    else # uploc == 'L'
1409
        if tfun === identity
1,194✔
1410
            for i in 1:m
452✔
1411
                for j in 1:n
3,836✔
1412
                    Cij = A[i,j] * (unit ? oB : B[j,j])
51,774✔
1413
                    for k in j + 1:n
37,840✔
1414
                        Cij += A[i,k] * B[k,j]
137,856✔
1415
                    end
245,544✔
1416
                    C[i,j] = Cij
34,004✔
1417
                end
64,172✔
1418
            end
3,836✔
1419
        else # tfun in (transpose, adjoint)
1420
            for i in 1:m
742✔
1421
                for j in n:-1:1
12,300✔
1422
                    Cij = A[i,j] * (unit ? oB : tfun(B[j,j]))
73,452✔
1423
                    for k in 1:j - 1
53,310✔
1424
                        Cij += A[i,k] * tfun(B[j,k])
211,140✔
1425
                    end
375,120✔
1426
                    C[i,j] = Cij
53,310✔
1427
                end
100,470✔
1428
            end
6,150✔
1429
        end
1430
    end
1431
    return C
2,273✔
1432
end
1433
# conjugate cases
1434
function generic_mattrimul!(C::AbstractMatrix, uploc, isunitc, ::Function, A::AbstractMatrix, xB::AdjOrTrans)
×
1435
    B = parent(xB)
×
1436
    require_one_based_indexing(C, A, B)
×
1437
    m, n = size(A, 1), size(A, 2)
×
1438
    N = size(B, 1)
×
1439
    if n != N
×
1440
        throw(DimensionMismatch(lazy"right hand side B needs first dimension of size $n, has size $N"))
×
1441
    end
1442
    mc, nc = size(C, 1), size(C, 2)
×
1443
    if mc != m || nc != N
×
1444
        throw(DimensionMismatch(lazy"output has dimensions ($mc,$nc), should have ($m,$N)"))
×
1445
    end
1446
    oB = oneunit(eltype(B))
×
1447
    unit = isunitc == 'U'
×
1448
    @inbounds if uploc == 'U'
×
1449
        for i in 1:m
×
1450
            for j in n:-1:1
×
1451
                Cij = A[i,j] * (unit ? oB : conj(B[j,j]))
×
1452
                for k in 1:j - 1
×
1453
                    Cij += A[i,k] * conj(B[k,j])
×
1454
                end
×
1455
                C[i,j] = Cij
×
1456
            end
×
1457
        end
×
1458
    else # uploc == 'L'
1459
        for i in 1:m
×
1460
            for j in 1:n
×
1461
                Cij = A[i,j] * (unit ? oB : conj(B[j,j]))
×
1462
                for k in j + 1:n
×
1463
                    Cij += A[i,k] * conj(B[k,j])
×
1464
                end
×
1465
                C[i,j] = Cij
×
1466
            end
×
1467
        end
×
1468
    end
1469
    return C
×
1470
end
1471

1472
#Generic solver using naive substitution
1473

1474
@inline _ustrip(a) = oneunit(a) \ a
322,882✔
1475
@inline _ustrip(a::Union{AbstractFloat,Integer,Complex,Rational}) = a
2,593,410✔
1476

1477
# manually hoisting b[j] significantly improves performance as of Dec 2015
1478
# manually eliding bounds checking significantly improves performance as of Dec 2015
1479
# replacing repeated references to A.data[j,j] with [Ajj = A.data[j,j] and references to Ajj]
1480
# does not significantly impact performance as of Dec 2015
1481
# in the transpose and conjugate transpose naive substitution variants,
1482
# accumulating in z rather than b[j,k] significantly improves performance as of Dec 2015
1483
function generic_trimatdiv!(C::AbstractVecOrMat, uploc, isunitc, tfun::Function, A::AbstractMatrix, B::AbstractVecOrMat)
19,326✔
1484
    require_one_based_indexing(C, A, B)
19,326✔
1485
    mA, nA = size(A)
19,326✔
1486
    m, n = size(B, 1), size(B,2)
19,326✔
1487
    if nA != m
19,326✔
1488
        throw(DimensionMismatch(lazy"second dimension of left hand side A, $nA, and first dimension of right hand side B, $m, must be equal"))
236✔
1489
    end
1490
    if size(C) != size(B)
38,180✔
1491
        throw(DimensionMismatch(lazy"size of output, $(size(C)), does not match size of right hand side, $(size(B))"))
×
1492
    end
1493
    iszero(mA) && return C
19,090✔
1494
    oA = oneunit(eltype(A))
19,084✔
1495
    @inbounds if uploc == 'U'
19,084✔
1496
        if isunitc == 'N'
10,182✔
1497
            if tfun === identity
8,847✔
1498
                for k in 1:n
6,056✔
1499
                    amm = A[m,m]
24,505✔
1500
                    iszero(amm) && throw(SingularException(m))
24,438✔
1501
                    Cm = C[m,k] = amm \ B[m,k]
24,358✔
1502
                    # fill C-column
1503
                    for i in m-1:-1:1
48,651✔
1504
                        C[i,k] = oA \ B[i,k] - _ustrip(A[i,m]) * Cm
202,709✔
1505
                    end
381,253✔
1506
                    for j in m-1:-1:1
48,651✔
1507
                        ajj = A[j,j]
203,306✔
1508
                        iszero(ajj) && throw(SingularException(j))
202,709✔
1509
                        Cj = C[j,k] = _ustrip(ajj) \ C[j,k]
202,709✔
1510
                        for i in j-1:-1:1
381,318✔
1511
                            C[i,k] -= _ustrip(A[i,j]) * Cj
845,369✔
1512
                        end
1,512,194✔
1513
                    end
381,253✔
1514
                end
42,740✔
1515
            else # tfun in (adjoint, transpose)
1516
                for k in 1:n
2,791✔
1517
                    for j in 1:m
14,635✔
1518
                        ajj = A[j,j]
136,654✔
1519
                        iszero(ajj) && throw(SingularException(j))
135,334✔
1520
                        Bj = B[j,k]
135,334✔
1521
                        for i in 1:j-1
135,334✔
1522
                            Bj -= tfun(A[i,j]) * C[i,k]
685,913✔
1523
                        end
1,002,643✔
1524
                        C[j,k] = tfun(ajj) \ Bj
135,334✔
1525
                    end
256,033✔
1526
                end
26,479✔
1527
            end
1528
        else # isunitc == 'U'
1529
            if tfun === identity
1,335✔
1530
                for k in 1:n
743✔
1531
                    Cm = C[m,k] = oA \ B[m,k]
5,044✔
1532
                    # fill C-column
1533
                    for i in m-1:-1:1
10,085✔
1534
                        C[i,k] = oA \ B[i,k] - _ustrip(A[i,m]) * Cm
41,671✔
1535
                    end
78,298✔
1536
                    for j in m-1:-1:1
10,085✔
1537
                        Cj = C[j,k]
41,671✔
1538
                        for i in 1:j-1
41,671✔
1539
                            C[i,k] -= _ustrip(A[i,j]) * Cj
151,868✔
1540
                        end
267,109✔
1541
                    end
78,298✔
1542
                end
5,044✔
1543
            else # tfun in (adjoint, transpose)
1544
                for k in 1:n
592✔
1545
                    for j in 1:m
1,840✔
1546
                        Bj = B[j,k]
16,920✔
1547
                        for i in 1:j-1
16,920✔
1548
                            Bj -= tfun(A[i,j]) * C[i,k]
87,956✔
1549
                        end
123,880✔
1550
                        C[j,k] = oA \ Bj
16,920✔
1551
                    end
32,000✔
1552
                end
1,840✔
1553
            end
1554
        end
1555
    else # uploc == 'L'
1556
        if isunitc == 'N'
8,902✔
1557
            if tfun === identity
7,571✔
1558
                for k in 1:n
5,131✔
1559
                    a11 = A[1,1]
20,844✔
1560
                    iszero(a11) && throw(SingularException(1))
20,778✔
1561
                    C1 = C[1,k] = a11 \ B[1,k]
20,778✔
1562
                    # fill C-column
1563
                    for i in 2:m
20,778✔
1564
                        C[i,k] = oA \ B[i,k] - _ustrip(A[i,1]) * C1
175,048✔
1565
                    end
329,318✔
1566
                    for j in 2:m
20,778✔
1567
                        ajj = A[j,j]
175,642✔
1568
                        iszero(ajj) && throw(SingularException(j))
175,048✔
1569
                        Cj = C[j,k] = _ustrip(ajj) \ C[j,k]
175,048✔
1570
                        for i in j+1:m
195,826✔
1571
                            C[i,k] -= _ustrip(A[i,j]) * Cj
742,730✔
1572
                        end
1,331,190✔
1573
                    end
329,318✔
1574
                end
20,778✔
1575
            else # tfun in (adjoint, transpose)
1576
                for k in 1:n
2,440✔
1577
                    for j in m:-1:1
26,062✔
1578
                        ajj = A[j,j]
121,374✔
1579
                        iszero(ajj) && throw(SingularException(j))
120,054✔
1580
                        Bj = B[j,k]
120,079✔
1581
                        for i in j+1:m
133,085✔
1582
                            Bj -= tfun(A[i,j]) * C[i,k]
617,689✔
1583
                        end
882,599✔
1584
                        C[j,k] = tfun(ajj) \ Bj
120,054✔
1585
                    end
227,077✔
1586
                end
13,031✔
1587
            end
1588
        else # isunitc == 'U'
1589
            if tfun === identity
1,331✔
1590
                for k in 1:n
739✔
1591
                    C1 = C[1,k] = oA \ B[1,k]
5,004✔
1592
                    # fill C-column
1593
                    for i in 2:m
5,004✔
1594
                        C[i,k] = oA \ B[i,k] - _ustrip(A[i,1]) * C1
41,311✔
1595
                    end
77,618✔
1596
                    for j in 2:m
5,004✔
1597
                        Cj = C[j,k]
41,311✔
1598
                        for i in j+1:m
46,315✔
1599
                            C[i,k] -= _ustrip(A[i,j]) * Cj
150,428✔
1600
                        end
264,549✔
1601
                    end
77,618✔
1602
                end
5,004✔
1603
            else # tfun in (adjoint, transpose)
1604
                for k in 1:n
592✔
1605
                    for j in m:-1:1
3,680✔
1606
                        Bj = B[j,k]
16,920✔
1607
                        for i in j+1:m
18,760✔
1608
                            Bj -= tfun(A[i,j]) * C[i,k]
87,956✔
1609
                        end
123,880✔
1610
                        C[j,k] = oA \ Bj
16,920✔
1611
                    end
32,000✔
1612
                end
1,840✔
1613
            end
1614
        end
1615
    end
1616
    return C
19,004✔
1617
end
1618
# conjugate cases
1619
function generic_trimatdiv!(C::AbstractVecOrMat, uploc, isunitc, ::Function, xA::AdjOrTrans, B::AbstractVecOrMat)
166✔
1620
    A = parent(xA)
166✔
1621
    require_one_based_indexing(C, A, B)
166✔
1622
    mA, nA = size(A)
166✔
1623
    m, n = size(B, 1), size(B,2)
166✔
1624
    if nA != m
166✔
1625
        throw(DimensionMismatch(lazy"second dimension of left hand side A, $nA, and first dimension of right hand side B, $m, must be equal"))
×
1626
    end
1627
    if size(C) != size(B)
332✔
1628
        throw(DimensionMismatch(lazy"size of output, $(size(C)), does not match size of right hand side, $(size(B))"))
×
1629
    end
1630
    iszero(mA) && return C
166✔
1631
    oA = oneunit(eltype(A))
164✔
1632
    @inbounds if uploc == 'U'
164✔
1633
        if isunitc == 'N'
82✔
1634
            for k in 1:n
82✔
1635
                amm = conj(A[m,m])
742✔
1636
                iszero(amm) && throw(SingularException(m))
742✔
1637
                Cm = C[m,k] = amm \ B[m,k]
742✔
1638
                # fill C-column
1639
                for i in m-1:-1:1
1,484✔
1640
                    C[i,k] = oA \ B[i,k] - _ustrip(conj(A[i,m])) * Cm
5,976✔
1641
                end
11,210✔
1642
                for j in m-1:-1:1
1,484✔
1643
                    ajj = conj(A[j,j])
5,976✔
1644
                    iszero(ajj) && throw(SingularException(j))
5,976✔
1645
                    Cj = C[j,k] = _ustrip(ajj) \ C[j,k]
5,976✔
1646
                    for i in j-1:-1:1
11,210✔
1647
                        C[i,k] -= _ustrip(conj(A[i,j])) * Cj
21,096✔
1648
                    end
36,958✔
1649
                end
11,210✔
1650
            end
1,402✔
1651
        else # isunitc == 'U'
1652
            for k in 1:n
×
1653
                Cm = C[m,k] = oA \ B[m,k]
×
1654
                # fill C-column
1655
                for i in m-1:-1:1
×
1656
                    C[i,k] = oA \ B[i,k] - _ustrip(conj(A[i,m])) * Cm
×
1657
                end
×
1658
                for j in m-1:-1:1
×
1659
                    Cj = C[j,k]
×
1660
                    for i in 1:j-1
×
1661
                        C[i,k] -= _ustrip(conj(A[i,j])) * Cj
×
1662
                    end
×
1663
                end
×
1664
            end
×
1665
        end
1666
    else # uploc == 'L'
1667
        if isunitc == 'N'
82✔
1668
            for k in 1:n
82✔
1669
                a11 = conj(A[1,1])
742✔
1670
                iszero(a11) && throw(SingularException(1))
742✔
1671
                C1 = C[1,k] = a11 \ B[1,k]
742✔
1672
                # fill C-column
1673
                for i in 2:m
742✔
1674
                    C[i,k] = oA \ B[i,k] - _ustrip(conj(A[i,1])) * C1
5,976✔
1675
                end
11,210✔
1676
                for j in 2:m
742✔
1677
                    ajj = conj(A[j,j])
5,976✔
1678
                    iszero(ajj) && throw(SingularException(j))
5,976✔
1679
                    Cj = C[j,k] = _ustrip(ajj) \ C[j,k]
5,976✔
1680
                    for i in j+1:m
6,718✔
1681
                        C[i,k] -= _ustrip(conj(A[i,j])) * Cj
21,096✔
1682
                    end
36,958✔
1683
                end
11,210✔
1684
            end
742✔
1685
        else # isunitc == 'U'
1686
            for k in 1:n
×
1687
                C1 = C[1,k] = oA \ B[1,k]
×
1688
                # fill C-column
1689
                for i in 2:m
×
1690
                    C[i,k] = oA \ B[i,k] - _ustrip(conj(A[i,1])) * C1
×
1691
                end
×
1692
                for j in 1:m
×
1693
                    Cj = C[j,k]
×
1694
                    for i in j+1:m
×
1695
                        C[i,k] -= _ustrip(conj(A[i,j])) * Cj
×
1696
                    end
×
1697
                end
×
1698
            end
×
1699
        end
1700
    end
1701
    return C
164✔
1702
end
1703

1704
function generic_mattridiv!(C::AbstractMatrix, uploc, isunitc, tfun::Function, A::AbstractMatrix, B::AbstractMatrix)
4,858✔
1705
    require_one_based_indexing(C, A, B)
4,858✔
1706
    m, n = size(A)
4,858✔
1707
    if size(B, 1) != n
4,858✔
1708
        throw(DimensionMismatch(lazy"right hand side B needs first dimension of size $n, has size $(size(B,1))"))
2,028✔
1709
    end
1710
    if size(C) != size(A)
5,660✔
1711
        throw(DimensionMismatch(lazy"size of output, $(size(C)), does not match size of left hand side, $(size(A))"))
×
1712
    end
1713
    oB = oneunit(eltype(B))
2,830✔
1714
    unit = isunitc == 'U'
2,830✔
1715
    @inbounds if uploc == 'U'
2,830✔
1716
        if tfun === identity
1,420✔
1717
            for i in 1:m
1,072✔
1718
                for j in 1:n
9,898✔
1719
                    Aij = A[i,j]
91,682✔
1720
                    for k in 1:j - 1
91,682✔
1721
                        Aij -= C[i,k]*B[k,j]
473,944✔
1722
                    end
678,112✔
1723
                    unit || (iszero(B[j,j]) && throw(SingularException(j)))
140,093✔
1724
                    C[i,j] = Aij / (unit ? oB : B[j,j])
91,682✔
1725
                end
173,466✔
1726
            end
18,724✔
1727
        else # tfun in (adjoint, transpose)
1728
            for i in 1:m
348✔
1729
                for j in n:-1:1
6,328✔
1730
                    Aij = A[i,j]
28,796✔
1731
                    for k in j + 1:n
31,960✔
1732
                        Aij -= C[i,k]*tfun(B[j,k])
142,992✔
1733
                    end
207,936✔
1734
                    unit || (iszero(B[j,j]) && throw(SingularException(j)))
44,956✔
1735
                    C[i,j] = Aij / (unit ? oB : tfun(B[j,j]))
28,796✔
1736
                end
54,428✔
1737
            end
5,980✔
1738
        end
1739
    else # uploc == 'L'
1740
        if tfun === identity
1,410✔
1741
            for i in 1:m
1,062✔
1742
                for j in n:-1:1
19,616✔
1743
                    Aij = A[i,j]
90,832✔
1744
                    for k in j + 1:n
100,640✔
1745
                        Aij -= C[i,k]*B[k,j]
470,244✔
1746
                    end
671,472✔
1747
                    unit || (iszero(B[j,j]) && throw(SingularException(j)))
138,793✔
1748
                    C[i,j] = Aij / (unit ? oB : B[j,j])
90,832✔
1749
                end
171,856✔
1750
            end
9,808✔
1751
        else # tfun in (adjoint, transpose)
1752
            for i in 1:m
348✔
1753
                for j in 1:n
3,164✔
1754
                    Aij = A[i,j]
28,796✔
1755
                    for k in 1:j - 1
28,796✔
1756
                        Aij -= C[i,k]*tfun(B[j,k])
142,992✔
1757
                    end
207,936✔
1758
                    unit || (iszero(B[j,j]) && throw(SingularException(j)))
44,956✔
1759
                    C[i,j] = Aij / (unit ? oB : tfun(B[j,j]))
28,796✔
1760
                end
54,428✔
1761
            end
3,164✔
1762
        end
1763
    end
1764
    return C
2,830✔
1765
end
1766
function generic_mattridiv!(C::AbstractMatrix, uploc, isunitc, ::Function, A::AbstractMatrix, xB::AdjOrTrans)
2✔
1767
    B = parent(xB)
2✔
1768
    require_one_based_indexing(C, A, B)
2✔
1769
    m, n = size(A)
2✔
1770
    if size(B, 1) != n
2✔
1771
        throw(DimensionMismatch(lazy"right hand side B needs first dimension of size $n, has size $(size(B,1))"))
×
1772
    end
1773
    if size(C) != size(A)
4✔
1774
        throw(DimensionMismatch(lazy"size of output, $(size(C)), does not match size of left hand side, $(size(A))"))
×
1775
    end
1776
    oB = oneunit(eltype(B))
2✔
1777
    unit = isunitc == 'U'
2✔
1778
    if uploc == 'U'
2✔
1779
        @inbounds for i in 1:m
2✔
1780
            for j in 1:n
20✔
1781
                Aij = A[i,j]
×
1782
                for k in 1:j - 1
×
1783
                    Aij -= C[i,k]*conj(B[k,j])
×
1784
                end
×
1785
                unit || (iszero(B[j,j]) && throw(SingularException(j)))
×
1786
                C[i,j] = Aij / (unit ? oB : conj(B[j,j]))
×
1787
            end
×
1788
        end
38✔
1789
    else # uploc == 'L'
1790
        @inbounds for i in 1:m
×
1791
            for j in n:-1:1
×
1792
                Aij = A[i,j]
×
1793
                for k in j + 1:n
×
1794
                    Aij -= C[i,k]*conj(B[k,j])
×
1795
                end
×
1796
                unit || (iszero(B[j,j]) && throw(SingularException(j)))
×
1797
                C[i,j] = Aij / (unit ? oB : conj(B[j,j]))
×
1798
            end
×
1799
        end
×
1800
    end
1801
    return C
2✔
1802
end
1803

1804
# these are needed because we don't keep track of left- and right-multiplication in tritrimul!
1805
rmul!(A::UpperTriangular, B::UpperTriangular)     = UpperTriangular(rmul!(triu!(A.data), B))
12✔
1806
rmul!(A::UpperTriangular, B::UnitUpperTriangular) = UpperTriangular(rmul!(triu!(A.data), B))
12✔
1807
rmul!(A::LowerTriangular, B::LowerTriangular)     = LowerTriangular(rmul!(tril!(A.data), B))
12✔
1808
rmul!(A::LowerTriangular, B::UnitLowerTriangular) = LowerTriangular(rmul!(tril!(A.data), B))
12✔
1809

1810
# Promotion
1811
## Promotion methods in matmul don't apply to triangular multiplication since
1812
## it is inplace. Hence we have to make very similar definitions, but without
1813
## allocation of a result array. For multiplication and unit diagonal division
1814
## the element type doesn't have to be stable under division whereas that is
1815
## necessary in the general triangular solve problem.
1816

1817
_inner_type_promotion(op, ::Type{TA}, ::Type{TB}) where {TA<:Integer,TB<:Integer} =
404✔
1818
    promote_op(matprod, TA, TB)
1819
_inner_type_promotion(op, ::Type{TA}, ::Type{TB}) where {TA,TB} =
6,368✔
1820
    promote_op(op, TA, TB)
1821
## The general promotion methods
1822
for mat in (:AbstractVector, :AbstractMatrix)
1823
    ### Left division with triangle to the left hence rhs cannot be transposed. No quotients.
1824
    @eval function \(A::Union{UnitUpperTriangular,UnitLowerTriangular}, B::$mat)
3,677✔
1825
        require_one_based_indexing(B)
4,235✔
1826
        TAB = _inner_type_promotion(\, eltype(A), eltype(B))
4,235✔
1827
        ldiv!(similar(B, TAB, size(B)), A, B)
4,793✔
1828
    end
1829
    ### Left division with triangle to the left hence rhs cannot be transposed. Quotients.
1830
    @eval function \(A::Union{UpperTriangular,LowerTriangular}, B::$mat)
4,350✔
1831
        require_one_based_indexing(B)
10,725✔
1832
        TAB = promote_op(\, eltype(A), eltype(B))
10,725✔
1833
        ldiv!(similar(B, TAB, size(B)), A, B)
17,100✔
1834
    end
1835
    ### Right division with triangle to the right hence lhs cannot be transposed. No quotients.
1836
    @eval function /(A::$mat, B::Union{UnitUpperTriangular, UnitLowerTriangular})
1,979✔
1837
        require_one_based_indexing(A)
2,537✔
1838
        TAB = _inner_type_promotion(/, eltype(A), eltype(B))
2,537✔
1839
        _rdiv!(similar(A, TAB, size(A)), A, B)
3,095✔
1840
    end
1841
    ### Right division with triangle to the right hence lhs cannot be transposed. Quotients.
1842
    @eval function /(A::$mat, B::Union{UpperTriangular,LowerTriangular})
2,005✔
1843
        require_one_based_indexing(A)
2,603✔
1844
        TAB = promote_op(/, eltype(A), eltype(B))
2,603✔
1845
        _rdiv!(similar(A, TAB, size(A)), A, B)
3,201✔
1846
    end
1847
end
1848

1849
## Some Triangular-Triangular cases. We might want to write tailored methods
1850
## for these cases, but I'm not sure it is worth it.
1851
for f in (:*, :\)
1852
    @eval begin
1853
        ($f)(A::LowerTriangular, B::LowerTriangular) =
743✔
1854
            LowerTriangular(@invoke $f(A::LowerTriangular, B::AbstractMatrix))
1855
        ($f)(A::LowerTriangular, B::UnitLowerTriangular) =
672✔
1856
            LowerTriangular(@invoke $f(A::LowerTriangular, B::AbstractMatrix))
1857
        ($f)(A::UnitLowerTriangular, B::LowerTriangular) =
627✔
1858
            LowerTriangular(@invoke $f(A::UnitLowerTriangular, B::AbstractMatrix))
1859
        ($f)(A::UnitLowerTriangular, B::UnitLowerTriangular) =
618✔
1860
            UnitLowerTriangular(@invoke $f(A::UnitLowerTriangular, B::AbstractMatrix))
1861
        ($f)(A::UpperTriangular, B::UpperTriangular) =
3,218✔
1862
            UpperTriangular(@invoke $f(A::UpperTriangular, B::AbstractMatrix))
1863
        ($f)(A::UpperTriangular, B::UnitUpperTriangular) =
672✔
1864
            UpperTriangular(@invoke $f(A::UpperTriangular, B::AbstractMatrix))
1865
        ($f)(A::UnitUpperTriangular, B::UpperTriangular) =
627✔
1866
            UpperTriangular(@invoke $f(A::UnitUpperTriangular, B::AbstractMatrix))
1867
        ($f)(A::UnitUpperTriangular, B::UnitUpperTriangular) =
623✔
1868
            UnitUpperTriangular(@invoke $f(A::UnitUpperTriangular, B::AbstractMatrix))
1869
    end
1870
end
1871
(/)(A::LowerTriangular, B::LowerTriangular) =
123✔
1872
    LowerTriangular(@invoke /(A::AbstractMatrix, B::LowerTriangular))
1873
(/)(A::LowerTriangular, B::UnitLowerTriangular) =
117✔
1874
    LowerTriangular(@invoke /(A::AbstractMatrix, B::UnitLowerTriangular))
1875
(/)(A::UnitLowerTriangular, B::LowerTriangular) =
98✔
1876
    LowerTriangular(@invoke /(A::AbstractMatrix, B::LowerTriangular))
1877
(/)(A::UnitLowerTriangular, B::UnitLowerTriangular) =
106✔
1878
    UnitLowerTriangular(@invoke /(A::AbstractMatrix, B::UnitLowerTriangular))
1879
(/)(A::UpperTriangular, B::UpperTriangular) =
167✔
1880
    UpperTriangular(@invoke /(A::AbstractMatrix, B::UpperTriangular))
1881
(/)(A::UpperTriangular, B::UnitUpperTriangular) =
117✔
1882
    UpperTriangular(@invoke /(A::AbstractMatrix, B::UnitUpperTriangular))
1883
(/)(A::UnitUpperTriangular, B::UpperTriangular) =
98✔
1884
    UpperTriangular(@invoke /(A::AbstractMatrix, B::UpperTriangular))
1885
(/)(A::UnitUpperTriangular, B::UnitUpperTriangular) =
106✔
1886
    UnitUpperTriangular(@invoke /(A::AbstractMatrix, B::UnitUpperTriangular))
1887

1888
# Complex matrix power for upper triangular factor, see:
1889
#   Higham and Lin, "A Schur-Padé algorithm for fractional powers of a Matrix",
1890
#     SIAM J. Matrix Anal. & Appl., 32 (3), (2011) 1056–1078.
1891
#   Higham and Lin, "An improved Schur-Padé algorithm for fractional powers of
1892
#     a matrix and their Fréchet derivatives", SIAM. J. Matrix Anal. & Appl.,
1893
#     34(3), (2013) 1341–1360.
1894
function powm!(A0::UpperTriangular, p::Real)
64✔
1895
    if abs(p) >= 1
64✔
1896
        throw(ArgumentError(lazy"p must be a real number in (-1,1), got $p"))
2✔
1897
    end
1898

1899
    normA0 = opnorm(A0, 1)
62✔
1900
    rmul!(A0, 1/normA0)
62✔
1901

1902
    theta = [1.53e-5, 2.25e-3, 1.92e-2, 6.08e-2, 1.25e-1, 2.03e-1, 2.84e-1]
434✔
1903
    n = checksquare(A0)
62✔
1904

1905
    A, m, s = invsquaring(A0, theta)
62✔
1906
    A = I - A
62✔
1907

1908
    # Compute accurate diagonal of I - T
1909
    sqrt_diag!(A0, A, s)
213✔
1910
    for i = 1:n
62✔
1911
        A[i, i] = -A[i, i]
214✔
1912
    end
366✔
1913
    # Compute the Padé approximant
1914
    c = 0.5 * (p - m) / (2 * m - 1)
62✔
1915
    triu!(A)
63✔
1916
    S = c * A
63✔
1917
    Stmp = similar(S)
62✔
1918
    for j = m-1:-1:1
124✔
1919
        j4 = 4 * j
266✔
1920
        c = (-p - j) / (j4 + 2)
266✔
1921
        for i = 1:n
266✔
1922
            @inbounds S[i, i] = S[i, i] + 1
940✔
1923
        end
1,614✔
1924
        copyto!(Stmp, S)
266✔
1925
        mul!(S, A, c)
270✔
1926
        ldiv!(Stmp, S)
270✔
1927

1928
        c = (p - j) / (j4 - 2)
266✔
1929
        for i = 1:n
266✔
1930
            @inbounds S[i, i] = S[i, i] + 1
940✔
1931
        end
1,614✔
1932
        copyto!(Stmp, S)
266✔
1933
        mul!(S, A, c)
270✔
1934
        ldiv!(Stmp, S)
270✔
1935
    end
470✔
1936
    for i = 1:n
62✔
1937
        S[i, i] = S[i, i] + 1
214✔
1938
    end
366✔
1939
    copyto!(Stmp, S)
62✔
1940
    mul!(S, A, -p)
63✔
1941
    ldiv!(Stmp, S)
63✔
1942
    for i = 1:n
62✔
1943
        @inbounds S[i, i] = S[i, i] + 1
214✔
1944
    end
366✔
1945

1946
    blockpower!(A0, S, p/(2^s))
62✔
1947
    for m = 1:s
62✔
1948
        mul!(Stmp.data, S, S)
248✔
1949
        copyto!(S, Stmp)
246✔
1950
        blockpower!(A0, S, p/(2^(s-m)))
246✔
1951
    end
430✔
1952
    rmul!(S, normA0^p)
63✔
1953
    return S
62✔
1954
end
1955
powm(A::LowerTriangular, p::Real) = copy(transpose(powm!(copy(transpose(A)), p::Real)))
2✔
1956

1957
# Complex matrix logarithm for the upper triangular factor, see:
1958
#   Al-Mohy and Higham, "Improved inverse scaling and squaring algorithms for
1959
#     the matrix logarithm", SIAM J. Sci. Comput., 34(4), (2012), pp. C153–C169.
1960
#   Al-Mohy, Higham and Relton, "Computing the Frechet derivative of the matrix
1961
#     logarithm and estimating the condition number", SIAM J. Sci. Comput.,
1962
#     35(4), (2013), C394–C410.
1963
#
1964
# Based on the code available at http://eprints.ma.man.ac.uk/1851/02/logm.zip,
1965
# Copyright (c) 2011, Awad H. Al-Mohy and Nicholas J. Higham
1966
# Julia version relicensed with permission from original authors
1967
log(A::UpperTriangular{T}) where {T<:BlasFloat} = log_quasitriu(A)
256✔
1968
log(A::UnitUpperTriangular{T}) where {T<:BlasFloat} = log_quasitriu(A)
26✔
1969
log(A::LowerTriangular) = copy(transpose(log(copy(transpose(A)))))
4✔
1970
log(A::UnitLowerTriangular) = copy(transpose(log(copy(transpose(A)))))
4✔
1971

1972
function log_quasitriu(A0::AbstractMatrix{T}) where T<:BlasFloat
295✔
1973
    # allocate real A if log(A) will be real and complex A otherwise
1974
    n = checksquare(A0)
295✔
1975
    if isreal(A0) && (!istriu(A0) || !any(x -> real(x) < zero(real(T)), diag(A0)))
1,235✔
1976
        A = T <: Complex ? real(A0) : copy(A0)
90✔
1977
    else
1978
        A = T <: Complex ? copy(A0) : complex(A0)
205✔
1979
    end
1980
    if A0 isa UnitUpperTriangular
295✔
1981
        A = UpperTriangular(parent(A))
26✔
1982
        @inbounds for i in 1:n
26✔
1983
            A[i,i] = 1
234✔
1984
        end
442✔
1985
    end
1986
    Y0 = _log_quasitriu!(A0, A)
385✔
1987
    # return complex result for complex input
1988
    Y = T <: Complex ? complex(Y0) : Y0
320✔
1989

1990
    if A0 isa UpperTriangular || A0 isa UnitUpperTriangular
295✔
1991
        return UpperTriangular(Y)
282✔
1992
    else
1993
        return Y
13✔
1994
    end
1995
end
1996
# type-stable implementation of log_quasitriu
1997
# A is a copy of A0 that is overwritten while computing the result. It has the same eltype
1998
# as the result.
1999
function _log_quasitriu!(A0, A)
295✔
2000
    # Find Padé degree m and s while replacing A with A^(1/2^s)
2001
    m, s = _find_params_log_quasitriu!(A)
295✔
2002

2003
    # Compute accurate superdiagonal of A
2004
    _pow_superdiag_quasitriu!(A, A0, 0.5^s)
583✔
2005

2006
    # Compute accurate block diagonal of A
2007
    _sqrt_pow_diag_quasitriu!(A, A0, s)
295✔
2008

2009
    # Get the Gauss-Legendre quadrature points and weights
2010
    R = zeros(Float64, m, m)
10,102✔
2011
    for i = 1:m - 1
295✔
2012
        R[i,i+1] = i / sqrt((2 * i)^2 - 1)
1,399✔
2013
        R[i+1,i] = R[i,i+1]
1,399✔
2014
    end
2,512✔
2015
    x,V = eigen(R)
590✔
2016
    w = Vector{Float64}(undef, m)
590✔
2017
    for i = 1:m
295✔
2018
        x[i] = (x[i] + 1) / 2
1,694✔
2019
        w[i] = V[1,i]^2
1,694✔
2020
    end
3,093✔
2021

2022
    # Compute the Padé approximation
2023
    t = eltype(A)
295✔
2024
    n = size(A, 1)
295✔
2025
    Y = zeros(t, n, n)
10,278✔
2026
    B = similar(A)
308✔
2027
    for k = 1:m
295✔
2028
        B .= t(x[k]) .* A
1,755✔
2029
        @inbounds for i in 1:n
1,694✔
2030
            B[i,i] += 1
8,412✔
2031
        end
15,130✔
2032
        Y .+= t(w[k]) .* rdiv_quasitriu!(A, B)
3,388✔
2033
    end
3,093✔
2034

2035
    # Scale back
2036
    lmul!(2.0^s, Y)
583✔
2037

2038
    # Compute accurate diagonal and superdiagonal of log(A)
2039
    _log_diag_quasitriu!(Y, A0)
295✔
2040

2041
    return Y
295✔
2042
end
2043

2044
# Auxiliary functions for matrix logarithm and matrix power
2045

2046
# Find Padé degree m and s while replacing A with A^(1/2^s)
2047
#   Al-Mohy and Higham, "Improved inverse scaling and squaring algorithms for
2048
#     the matrix logarithm", SIAM J. Sci. Comput., 34(4), (2012), pp. C153–C169.
2049
#   from Algorithm 4.1
2050
function _find_params_log_quasitriu!(A)
295✔
2051
    maxsqrt = 100
295✔
2052
    theta = [1.586970738772063e-005,
2,065✔
2053
         2.313807884242979e-003,
2054
         1.938179313533253e-002,
2055
         6.209171588994762e-002,
2056
         1.276404810806775e-001,
2057
         2.060962623452836e-001,
2058
         2.879093714241194e-001]
2059
    tmax = size(theta, 1)
295✔
2060
    n = size(A, 1)
295✔
2061
    p = 0
295✔
2062
    m = 0
295✔
2063

2064
    # Find s0, the smallest s such that the ρ(triu(A)^(1/2^s) - I) ≤ theta[tmax], where ρ(X)
2065
    # is the spectral radius of X
2066
    d = complex.(@view(A[diagind(A)]))
295✔
2067
    dm1 = d .- 1
590✔
2068
    s = 0
295✔
2069
    while norm(dm1, Inf) > theta[tmax] && s < maxsqrt
1,384✔
2070
        d .= sqrt.(d)
1,089✔
2071
        dm1 .= d .- 1
2,178✔
2072
        s = s + 1
1,089✔
2073
    end
1,089✔
2074
    s0 = s
295✔
2075

2076
    # Compute repeated roots
2077
    for k = 1:min(s, maxsqrt)
295✔
2078
        _sqrt_quasitriu!(A isa UpperTriangular ? parent(A) : A, A)
1,089✔
2079
    end
1,917✔
2080

2081
    # these three never needed at the same time, so reuse the same temporary
2082
    AmI = AmI4 = AmI5 = A - I
1,401✔
2083
    AmI2 = AmI * AmI
308✔
2084
    AmI3 = AmI2 * AmI
308✔
2085
    d2 = sqrt(opnorm(AmI2, 1))
295✔
2086
    d3 = cbrt(opnorm(AmI3, 1))
295✔
2087
    alpha2 = max(d2, d3)
295✔
2088
    foundm = false
295✔
2089
    if alpha2 <= theta[2]
295✔
2090
        m = alpha2 <= theta[1] ? 1 : 2
9✔
2091
        foundm = true
9✔
2092
    end
2093

2094
    while !foundm
633✔
2095
        more_sqrt = false
624✔
2096
        mul!(AmI4, AmI2, AmI2)
624✔
2097
        d4 = opnorm(AmI4, 1)^(1/4)
1,248✔
2098
        alpha3 = max(d3, d4)
624✔
2099
        if alpha3 <= theta[tmax]
624✔
2100
            local j
2101
            for outer j = 3:tmax
359✔
2102
                if alpha3 <= theta[j]
1,478✔
2103
                    break
359✔
2104
                end
2105
            end
1,119✔
2106
            if j <= 6
359✔
2107
                m = j
220✔
2108
                break
220✔
2109
            elseif alpha3 / 2 <= theta[5] && p < 2
139✔
2110
                more_sqrt = true
98✔
2111
                p = p + 1
98✔
2112
           end
2113
        end
2114

2115
        if !more_sqrt
404✔
2116
            mul!(AmI5, AmI3, AmI2)
306✔
2117
            d5 = opnorm(AmI5, 1)^(1/5)
612✔
2118
            alpha4 = max(d4, d5)
306✔
2119
            eta = min(alpha3, alpha4)
306✔
2120
            if eta <= theta[tmax]
306✔
2121
                j = 0
65✔
2122
                for outer j = 6:tmax
65✔
2123
                    if eta <= theta[j]
130✔
2124
                        m = j
65✔
2125
                        break
65✔
2126
                    end
2127
                end
65✔
2128
                break
65✔
2129
            end
2130
        end
2131

2132
        if s == maxsqrt
339✔
2133
            m = tmax
1✔
2134
            break
1✔
2135
        end
2136
        _sqrt_quasitriu!(A isa UpperTriangular ? parent(A) : A, A)
338✔
2137
        copyto!(AmI, A)
440✔
2138
        for i in 1:n
338✔
2139
            @inbounds AmI[i,i] -= 1
1,726✔
2140
        end
3,114✔
2141
        mul!(AmI2, AmI, AmI)
338✔
2142
        mul!(AmI3, AmI2, AmI)
338✔
2143
        d3 = cbrt(opnorm(AmI3, 1))
338✔
2144
        s = s + 1
338✔
2145
    end
338✔
2146
    return m, s
295✔
2147
end
2148

2149
# Compute accurate diagonal of A = A0^s - I
2150
function sqrt_diag!(A0::UpperTriangular, A::UpperTriangular, s)
1✔
2151
    n = checksquare(A0)
62✔
2152
    T = eltype(A)
62✔
2153
    @inbounds for i = 1:n
62✔
2154
        a = complex(A0[i,i])
214✔
2155
        A[i,i] = _sqrt_pow(a, s)
214✔
2156
    end
365✔
2157
end
2158
# Compute accurate block diagonal of A = A0^s - I for upper quasi-triangular A0 produced
2159
# by the Schur decomposition. Diagonal is made of 1x1 and 2x2 blocks.
2160
# 2x2 blocks are real with non-negative conjugate pair eigenvalues
2161
function _sqrt_pow_diag_quasitriu!(A, A0, s)
295✔
2162
    n = checksquare(A0)
295✔
2163
    t = typeof(sqrt(zero(eltype(A))))
295✔
2164
    i = 1
295✔
2165
    @inbounds while i < n
295✔
2166
        if iszero(A0[i+1,i])  # 1x1 block
1,131✔
2167
            A[i,i] = _sqrt_pow(t(A0[i,i]), s)
1,114✔
2168
            i += 1
1,114✔
2169
        else  # real 2x2 block
2170
            @views _sqrt_pow_diag_block_2x2!(A[i:i+1,i:i+1], A0[i:i+1,i:i+1], s)
17✔
2171
            i += 2
17✔
2172
        end
2173
    end
1,131✔
2174
    if i == n  # last block is 1x1
295✔
2175
        @inbounds A[n,n] = _sqrt_pow(t(A0[n,n]), s)
288✔
2176
    end
2177
    return A
295✔
2178
end
2179
# compute a^(1/2^s)-1
2180
#   Al-Mohy, "A more accurate Briggs method for the logarithm",
2181
#      Numer. Algorithms, 59, (2012), 393–402.
2182
#   Algorithm 2
2183
function _sqrt_pow(a::Number, s)
1,616✔
2184
    T = typeof(sqrt(zero(a)))
1,616✔
2185
    s == 0 && return T(a) - 1
1,616✔
2186
    s0 = s
1,589✔
2187
    if imag(a) >= 0 && real(a) <= 0 && !iszero(a)  # angle(a) ≥ π / 2
1,589✔
2188
        a = sqrt(a)
139✔
2189
        s0 = s - 1
139✔
2190
    end
2191
    z0 = a - 1
1,589✔
2192
    a = sqrt(a)
1,589✔
2193
    r = 1 + a
1,589✔
2194
    for j = 1:s0-1
1,589✔
2195
        a = sqrt(a)
4,444✔
2196
        r = r * (1 + a)
4,444✔
2197
    end
7,315✔
2198
    return z0 / r
1,589✔
2199
end
2200
# compute A0 = A^(1/2^s)-I for 2x2 real matrices A and A0
2201
# A has non-negative conjugate pair eigenvalues
2202
# "Improved Inverse Scaling and Squaring Algorithms for the Matrix Logarithm"
2203
# SIAM J. Sci. Comput., 34(4), (2012) C153–C169. doi: 10.1137/110852553
2204
# Algorithm 5.1
2205
Base.@propagate_inbounds function _sqrt_pow_diag_block_2x2!(A, A0, s)
2206
    _sqrt_real_2x2!(A, A0)
17✔
2207
    if isone(s)
17✔
2208
        A[1,1] -= 1
×
2209
        A[2,2] -= 1
×
2210
    else
2211
        # Z = A - I
2212
        z11, z21, z12, z22 = A[1,1] - 1, A[2,1], A[1,2], A[2,2] - 1
17✔
2213
        # A = sqrt(A)
2214
        _sqrt_real_2x2!(A, A)
17✔
2215
        # P = A + I
2216
        p11, p21, p12, p22 = A[1,1] + 1, A[2,1], A[1,2], A[2,2] + 1
17✔
2217
        for i in 1:(s - 2)
17✔
2218
            # A = sqrt(A)
2219
            _sqrt_real_2x2!(A, A)
423✔
2220
            a11, a21, a12, a22 = A[1,1], A[2,1], A[1,2], A[2,2]
423✔
2221
            # P += P * A
2222
            r11 = p11*(1 + a11) + p12*a21
423✔
2223
            r22 = p21*a12 + p22*(1 + a22)
423✔
2224
            p21 = p21*(1 + a11) + p22*a21
423✔
2225
            p12 = p11*a12 + p12*(1 + a22)
423✔
2226
            p11 = r11
423✔
2227
            p22 = r22
423✔
2228
        end
829✔
2229
        # A = Z / P
2230
        c = inv(p11*p22 - p21*p12)
17✔
2231
        A[1,1] = (p22*z11 - p21*z12) * c
17✔
2232
        A[2,1] = (p22*z21 - p21*z22) * c
17✔
2233
        A[1,2] = (p11*z12 - p12*z11) * c
17✔
2234
        A[2,2] = (p11*z22 - p12*z21) * c
17✔
2235
    end
2236
    return A
17✔
2237
end
2238
# Compute accurate superdiagonal of A = A0^s - I for upper quasi-triangular A0 produced
2239
# by a Schur decomposition.
2240
# Higham and Lin, "A Schur–Padé Algorithm for Fractional Powers of a Matrix"
2241
# SIAM J. Matrix Anal. Appl., 32(3), (2011), 1056–1078.
2242
# Equation 5.6
2243
# see also blockpower for when A0 is upper triangular
2244
function _pow_superdiag_quasitriu!(A, A0, p)
295✔
2245
    n = checksquare(A0)
295✔
2246
    t = eltype(A)
295✔
2247
    k = 1
295✔
2248
    @inbounds while k < n
295✔
2249
        if !iszero(A[k+1,k])
1,124✔
2250
            k += 2
10✔
2251
            continue
10✔
2252
        end
2253
        if !(k == n - 1 || iszero(A[k+2,k+1]))
1,949✔
2254
            k += 3
7✔
2255
            continue
7✔
2256
        end
2257
        Ak = t(A0[k,k])
1,107✔
2258
        Akp1 = t(A0[k+1,k+1])
1,107✔
2259

2260
        Akp = Ak^p
1,107✔
2261
        Akp1p = Akp1^p
1,107✔
2262

2263
        if Ak == Akp1
1,107✔
2264
            A[k,k+1] = p * A0[k,k+1] * Ak^(p-1)
306✔
2265
        elseif 2 * abs(Ak) < abs(Akp1) || 2 * abs(Akp1) < abs(Ak) || iszero(Akp1 + Ak)
1,590✔
2266
            A[k,k+1] = A0[k,k+1] * (Akp1p - Akp) / (Akp1 - Ak)
424✔
2267
        else
2268
            logAk = log(Ak)
377✔
2269
            logAkp1 = log(Akp1)
377✔
2270
            z = (Akp1 - Ak)/(Akp1 + Ak)
377✔
2271
            if abs(z) > 1
488✔
2272
                A[k,k+1] = A0[k,k+1] * (Akp1p - Akp) / (Akp1 - Ak)
73✔
2273
            else
2274
                w = atanh(z) + im * pi * (unw(logAkp1-logAk) - unw(log1p(z)-log1p(-z)))
399✔
2275
                dd = 2 * exp(p*(logAk+logAkp1)/2) * sinh(p*w) / (Akp1 - Ak);
513✔
2276
                A[k,k+1] = A0[k,k+1] * dd
304✔
2277
            end
2278
        end
2279
        k += 1
1,107✔
2280
    end
1,124✔
2281
end
2282

2283
# Compute accurate block diagonal and superdiagonal of A = log(A0) for upper
2284
# quasi-triangular A0 produced by the Schur decomposition.
2285
function _log_diag_quasitriu!(A, A0)
295✔
2286
    n = checksquare(A0)
295✔
2287
    t = eltype(A)
295✔
2288
    k = 1
295✔
2289
    @inbounds while k < n
295✔
2290
        if iszero(A0[k+1,k])  # 1x1 block
776✔
2291
            Ak = t(A0[k,k])
759✔
2292
            logAk = log(Ak)
759✔
2293
            A[k,k] = logAk
759✔
2294
            if k < n - 2 && iszero(A0[k+2,k+1])
759✔
2295
                Akp1 = t(A0[k+1,k+1])
355✔
2296
                logAkp1 = log(Akp1)
355✔
2297
                A[k+1,k+1] = logAkp1
355✔
2298
                if Ak == Akp1
355✔
2299
                    A[k,k+1] = A0[k,k+1] / Ak
114✔
2300
                elseif 2 * abs(Ak) < abs(Akp1) || 2 * abs(Akp1) < abs(Ak) || iszero(Akp1 + Ak)
486✔
2301
                    A[k,k+1] = A0[k,k+1] * (logAkp1 - logAk) / (Akp1 - Ak)
121✔
2302
                else
2303
                    z = (Akp1 - Ak)/(Akp1 + Ak)
120✔
2304
                    if abs(z) > 1
152✔
2305
                        A[k,k+1] = A0[k,k+1] * (logAkp1 - logAk) / (Akp1 - Ak)
32✔
2306
                    else
2307
                        w = atanh(z) + im * pi * (unw(logAkp1-logAk) - unw(log1p(z)-log1p(-z)))
128✔
2308
                        A[k,k+1] = 2 * A0[k,k+1] * w / (Akp1 - Ak)
88✔
2309
                    end
2310
                end
2311
                k += 2
355✔
2312
            else
2313
                k += 1
404✔
2314
            end
2315
        else  # real 2x2 block
2316
            @views _log_diag_block_2x2!(A[k:k+1,k:k+1], A0[k:k+1,k:k+1])
17✔
2317
            k += 2
17✔
2318
        end
2319
    end
776✔
2320
    if k == n  # last 1x1 block
295✔
2321
        @inbounds A[n,n] = log(t(A0[n,n]))
288✔
2322
    end
2323
    return A
295✔
2324
end
2325
# compute A0 = log(A) for 2x2 real matrices A and A0, where A0 is a diagonal 2x2 block
2326
# produced by real Schur decomposition.
2327
# Al-Mohy, Higham and Relton, "Computing the Frechet derivative of the matrix
2328
# logarithm and estimating the condition number", SIAM J. Sci. Comput.,
2329
# 35(4), (2013), C394–C410.
2330
# Eq. 6.1
2331
Base.@propagate_inbounds function _log_diag_block_2x2!(A, A0)
2332
    a, b, c = A0[1,1], A0[1,2], A0[2,1]
17✔
2333
    # avoid underflow/overflow for large/small b and c
2334
    s = sqrt(abs(b)) * sqrt(abs(c))
17✔
2335
    θ = atan(s, a)
17✔
2336
    t = θ / s
17✔
2337
    au = abs(a)
17✔
2338
    if au > s
17✔
2339
        a1 = log1p((s / au)^2) / 2 + log(au)
7✔
2340
    else
2341
        a1 = log1p((au / s)^2) / 2 + log(s)
10✔
2342
    end
2343
    A[1,1] = a1
17✔
2344
    A[2,1] = c*t
17✔
2345
    A[1,2] = b*t
17✔
2346
    A[2,2] = a1
17✔
2347
    return A
17✔
2348
end
2349

2350
# Used only by powm at the moment
2351
# Repeatedly compute the square roots of A so that in the end its
2352
# eigenvalues are close enough to the positive real line
2353
function invsquaring(A0::UpperTriangular, theta)
62✔
2354
    require_one_based_indexing(theta)
62✔
2355
    # assumes theta is in ascending order
2356
    maxsqrt = 100
62✔
2357
    tmax = size(theta, 1)
62✔
2358
    n = checksquare(A0)
62✔
2359
    A = complex(copy(A0))
62✔
2360
    p = 0
62✔
2361
    m = 0
62✔
2362

2363
    # Compute repeated roots
2364
    d = complex(diag(A))
62✔
2365
    dm1 = d .- 1
124✔
2366
    s = 0
62✔
2367
    while norm(dm1, Inf) > theta[tmax] && s < maxsqrt
216✔
2368
        d .= sqrt.(d)
154✔
2369
        dm1 .= d .- 1
308✔
2370
        s = s + 1
154✔
2371
    end
154✔
2372
    s0 = s
62✔
2373
    for k = 1:min(s, maxsqrt)
62✔
2374
        A = sqrt(A)
154✔
2375
    end
262✔
2376

2377
    AmI = A - I
213✔
2378
    d2 = sqrt(opnorm(AmI^2, 1))
62✔
2379
    d3 = cbrt(opnorm(AmI^3, 1))
62✔
2380
    alpha2 = max(d2, d3)
63✔
2381
    foundm = false
62✔
2382
    if alpha2 <= theta[2]
62✔
2383
        m = alpha2 <= theta[1] ? 1 : 2
×
2384
        foundm = true
×
2385
    end
2386

2387
    while !foundm
198✔
2388
        more = false
198✔
2389
        if s > s0
198✔
2390
            d3 = cbrt(opnorm(AmI^3, 1))
122✔
2391
        end
2392
        d4 = opnorm(AmI^4, 1)^(1/4)
392✔
2393
        alpha3 = max(d3, d4)
202✔
2394
        if alpha3 <= theta[tmax]
198✔
2395
            local j
2396
            for outer j = 3:tmax
150✔
2397
                if alpha3 <= theta[j]
644✔
2398
                    break
150✔
2399
                elseif alpha3 / 2 <= theta[5] && p < 2
494✔
2400
                    more = true
124✔
2401
                    p = p + 1
124✔
2402
                end
2403
            end
494✔
2404
            if j <= 6
150✔
2405
                m = j
62✔
2406
                foundm = true
62✔
2407
                break
62✔
2408
            elseif alpha3 / 2 <= theta[5] && p < 2
88✔
2409
                more = true
×
2410
                p = p + 1
×
2411
           end
2412
        end
2413

2414
        if !more
136✔
2415
            d5 = opnorm(AmI^5, 1)^(1/5)
182✔
2416
            alpha4 = max(d4, d5)
94✔
2417
            eta = min(alpha3, alpha4)
94✔
2418
            if eta <= theta[tmax]
92✔
2419
                j = 0
60✔
2420
                for outer j = 6:tmax
60✔
2421
                    if eta <= theta[j]
60✔
2422
                        m = j
14✔
2423
                        break
14✔
2424
                    end
2425
                    break
46✔
2426
                end
×
2427
            end
2428
            if s == maxsqrt
92✔
2429
                m = tmax
×
2430
                break
×
2431
            end
2432
            A = sqrt(A)
92✔
2433
            AmI = A - I
300✔
2434
            s = s + 1
92✔
2435
        end
2436
    end
136✔
2437

2438
    # Compute accurate superdiagonal of T
2439
    p = 1 / 2^s
62✔
2440
    A = complex(A)
62✔
2441
    blockpower!(A, A0, p)
62✔
2442
    return A,m,s
62✔
2443
end
2444

2445
# Compute accurate diagonal and superdiagonal of A = A0^p
2446
function blockpower!(A::UpperTriangular, A0::UpperTriangular, p)
370✔
2447
    n = checksquare(A0)
370✔
2448
    @inbounds for k = 1:n-1
370✔
2449
        Ak = complex(A0[k,k])
968✔
2450
        Akp1 = complex(A0[k+1,k+1])
968✔
2451

2452
        Akp = Ak^p
968✔
2453
        Akp1p = Akp1^p
968✔
2454

2455
        A[k,k] = Akp
968✔
2456
        A[k+1,k+1] = Akp1p
968✔
2457

2458
        if Ak == Akp1
968✔
2459
            A[k,k+1] = p * A0[k,k+1] * Ak^(p-1)
137✔
2460
        elseif 2 * abs(Ak) < abs(Akp1) || 2 * abs(Akp1) < abs(Ak) || iszero(Akp1 + Ak)
1,618✔
2461
            A[k,k+1] = A0[k,k+1] * (Akp1p - Akp) / (Akp1 - Ak)
102✔
2462
        else
2463
            logAk = log(Ak)
729✔
2464
            logAkp1 = log(Akp1)
729✔
2465
            z = (Akp1 - Ak)/(Akp1 + Ak)
729✔
2466
            if abs(z) > 1
729✔
2467
                A[k,k+1] = A0[k,k+1] * (Akp1p - Akp) / (Akp1 - Ak)
×
2468
            else
2469
                w = atanh(z) + im * pi * (unw(logAkp1-logAk) - unw(log1p(z)-log1p(-z)))
729✔
2470
                dd = 2 * exp(p*(logAk+logAkp1)/2) * sinh(p*w) / (Akp1 - Ak);
1,458✔
2471
                A[k,k+1] = A0[k,k+1] * dd
729✔
2472
            end
2473
        end
2474
    end
968✔
2475
end
2476

2477
# Unwinding number
2478
unw(x::Real) = 0
×
2479
unw(x::Number) = ceil((imag(x) - pi) / (2 * pi))
1,972✔
2480

2481
# compute A / B for upper quasi-triangular B, possibly overwriting B
2482
function rdiv_quasitriu!(A, B)
1,694✔
2483
    n = checksquare(A)
1,694✔
2484
    AG = copy(A)
1,694✔
2485
    # use Givens rotations to annihilate 2x2 blocks
2486
    @inbounds for k in 1:(n-1)
1,694✔
2487
        s = B[k+1,k]
6,718✔
2488
        iszero(s) && continue  # 1x1 block
6,718✔
2489
        G = first(givens(B[k+1,k+1], s, k, k+1))
85✔
2490
        rmul!(B, G)
400✔
2491
        rmul!(AG, G)
85✔
2492
    end
11,758✔
2493
    return rdiv!(AG, UpperTriangular(B))
1,694✔
2494
end
2495

2496
# End of auxiliary functions for matrix logarithm and matrix power
2497

2498
sqrt(A::UpperTriangular) = sqrt_quasitriu(A)
590✔
2499
function sqrt(A::UnitUpperTriangular{T}) where T
20✔
2500
    B = A.data
20✔
2501
    n = checksquare(B)
20✔
2502
    t = typeof(sqrt(zero(T)))
20✔
2503
    R = Matrix{t}(I, n, n)
20✔
2504
    tt = typeof(oneunit(t)*oneunit(t))
20✔
2505
    half = inv(R[1,1]+R[1,1]) # for general, algebraic cases. PR#20214
22✔
2506
    @inbounds for j = 1:n
20✔
2507
        for i = j-1:-1:1
264✔
2508
            r::tt = B[i,j]
520✔
2509
            @simd for k = i+1:j-1
640✔
2510
                r -= R[i,k]*R[k,j]
1,180✔
2511
            end
2512
            iszero(r) || (R[i,j] = half*r)
1,035✔
2513
        end
914✔
2514
    end
264✔
2515
    return UnitUpperTriangular(R)
20✔
2516
end
2517
sqrt(A::LowerTriangular) = copy(transpose(sqrt(copy(transpose(A)))))
8✔
2518
sqrt(A::UnitLowerTriangular) = copy(transpose(sqrt(copy(transpose(A)))))
8✔
2519

2520
# Auxiliary functions for matrix square root
2521

2522
# square root of upper triangular or real upper quasitriangular matrix
2523
function sqrt_quasitriu(A0; blockwidth = eltype(A0) <: Complex ? 512 : 256)
1,356✔
2524
    n = checksquare(A0)
678✔
2525
    T = eltype(A0)
678✔
2526
    Tr = typeof(sqrt(real(zero(T))))
678✔
2527
    Tc = typeof(sqrt(complex(zero(T))))
678✔
2528
    if isreal(A0)
1,981✔
2529
        is_sqrt_real = true
288✔
2530
        if istriu(A0)
288✔
2531
            for i in 1:n
204✔
2532
                Aii = real(A0[i,i])
592✔
2533
                if Aii < zero(Aii)
592✔
2534
                    is_sqrt_real = false
40✔
2535
                    break
40✔
2536
                end
2537
            end
552✔
2538
        end
2539
        if is_sqrt_real
288✔
2540
            R = zeros(Tr, n, n)
26,860✔
2541
            A = real(A0)
248✔
2542
        else
2543
            R = zeros(Tc, n, n)
347✔
2544
            A = A0
40✔
2545
        end
2546
    else
2547
        A = A0
390✔
2548
        R = zeros(Tc, n, n)
390✔
2549
    end
2550
    _sqrt_quasitriu!(R, A; blockwidth=blockwidth, n=n)
921✔
2551
    Rc = eltype(A0) <: Real ? R : complex(R)
775✔
2552
    if A0 isa UpperTriangular
678✔
2553
        return UpperTriangular(Rc)
590✔
2554
    elseif A0 isa UnitUpperTriangular
88✔
2555
        return UnitUpperTriangular(Rc)
×
2556
    else
2557
        return Rc
88✔
2558
    end
2559
end
2560

2561
# in-place recursive sqrt of upper quasi-triangular matrix A from
2562
# Deadman E., Higham N.J., Ralha R. (2013) Blocked Schur Algorithms for Computing the Matrix
2563
# Square Root. Applied Parallel and Scientific Computing. PARA 2012. Lecture Notes in
2564
# Computer Science, vol 7782. https://doi.org/10.1007/978-3-642-36803-5_12
2565
function _sqrt_quasitriu!(R, A; blockwidth=64, n=checksquare(A))
4,326✔
2566
    if n ≤ blockwidth || !(eltype(R) <: BlasFloat) # base case, perform "point" algorithm
2,199✔
2567
        _sqrt_quasitriu_block!(R, A)
2,136✔
2568
    else  # compute blockwise recursion
2569
        split = div(n, 2)
31✔
2570
        iszero(A[split+1, split]) || (split += 1) # don't split 2x2 diagonal block
39✔
2571
        r1 = 1:split
31✔
2572
        r2 = (split + 1):n
31✔
2573
        n1, n2 = split, n - split
31✔
2574
        A11, A12, A22 = @views A[r1,r1], A[r1,r2], A[r2,r2]
31✔
2575
        R11, R12, R22 = @views R[r1,r1], R[r1,r2], R[r2,r2]
31✔
2576
        # solve diagonal blocks recursively
2577
        _sqrt_quasitriu!(R11, A11; blockwidth=blockwidth, n=n1)
31✔
2578
        _sqrt_quasitriu!(R22, A22; blockwidth=blockwidth, n=n2)
31✔
2579
        # solve off-diagonal block
2580
        R12 .= .- A12
62✔
2581
        _sylvester_quasitriu!(R11, R22, R12; blockwidth=blockwidth, nA=n1, nB=n2, raise=false)
31✔
2582
    end
2583
    return R
2,167✔
2584
end
2585

2586
function _sqrt_quasitriu_block!(R, A)
2587
    _sqrt_quasitriu_diag_block!(R, A)
2,136✔
2588
    _sqrt_quasitriu_offdiag_block!(R, A)
2,136✔
2589
    return R
2,136✔
2590
end
2591

2592
function _sqrt_quasitriu_diag_block!(R, A)
2,136✔
2593
    n = size(R, 1)
2,136✔
2594
    ta = eltype(R) <: Complex ? complex(eltype(A)) : eltype(A)
2,136✔
2595
    i = 1
2,136✔
2596
    @inbounds while i < n
2,136✔
2597
        if iszero(A[i + 1, i])
7,770✔
2598
            R[i, i] = sqrt(ta(A[i, i]))
7,034✔
2599
            i += 1
7,034✔
2600
        else
2601
            # This branch is never reached when A is complex triangular
2602
            @assert eltype(A) <: Real
736✔
2603
            @views _sqrt_real_2x2!(R[i:(i + 1), i:(i + 1)], A[i:(i + 1), i:(i + 1)])
736✔
2604
            i += 2
736✔
2605
        end
2606
    end
7,770✔
2607
    if i == n
2,136✔
2608
        R[n, n] = sqrt(ta(A[n, n]))
1,698✔
2609
    end
2610
    return R
2,136✔
2611
end
2612

2613
function _sqrt_quasitriu_offdiag_block!(R, A)
2,136✔
2614
    n = size(R, 1)
2,136✔
2615
    j = 1
2,136✔
2616
    @inbounds while j ≤ n
2,136✔
2617
        jsize_is_2 = j < n && !iszero(A[j + 1, j])
9,468✔
2618
        i = j - 1
9,468✔
2619
        while i > 0
76,283✔
2620
            isize_is_2 = i > 1 && !iszero(A[i, i - 1])
66,815✔
2621
            if isize_is_2
66,815✔
2622
                if jsize_is_2
831✔
2623
                    _sqrt_quasitriu_offdiag_block_2x2!(R, A, i - 1, j)
417✔
2624
                else
2625
                    _sqrt_quasitriu_offdiag_block_2x1!(R, A, i - 1, j)
417✔
2626
                end
2627
                i -= 2
831✔
2628
            else
2629
                if jsize_is_2
65,984✔
2630
                    _sqrt_quasitriu_offdiag_block_1x2!(R, A, i, j)
250✔
2631
                else
2632
                    _sqrt_quasitriu_offdiag_block_1x1!(R, A, i, j)
65,734✔
2633
                end
2634
                i -= 1
65,984✔
2635
            end
2636
        end
66,815✔
2637
        j += 2 - !jsize_is_2
9,468✔
2638
    end
9,468✔
2639
    return R
2,136✔
2640
end
2641

2642
# real square root of 2x2 diagonal block of quasi-triangular matrix from real Schur
2643
# decomposition. Eqs 6.8-6.9 and Algorithm 6.5 of
2644
# Higham, 2008, "Functions of Matrices: Theory and Computation", SIAM.
2645
Base.@propagate_inbounds function _sqrt_real_2x2!(R, A)
2646
    # in the real Schur form, A[1, 1] == A[2, 2], and A[2, 1] * A[1, 2] < 0
2647
    θ, a21, a12 = A[1, 1], A[2, 1], A[1, 2]
1,193✔
2648
    # avoid overflow/underflow of μ
2649
    # for real sqrt, |d| ≤ 2 max(|a12|,|a21|)
2650
    μ = sqrt(abs(a12)) * sqrt(abs(a21))
1,193✔
2651
    α = _real_sqrt(θ, μ)
1,307✔
2652
    c = 2α
1,193✔
2653
    R[1, 1] = α
1,193✔
2654
    R[2, 1] = a21 / c
1,193✔
2655
    R[1, 2] = a12 / c
1,193✔
2656
    R[2, 2] = α
1,193✔
2657
    return R
1,193✔
2658
end
2659

2660
# real part of square root of θ+im*μ
2661
@inline function _real_sqrt(θ, μ)
2662
    t = sqrt((abs(θ) + hypot(θ, μ)) / 2)
1,278✔
2663
    return θ ≥ 0 ? t : μ / 2t
1,307✔
2664
end
2665

2666
Base.@propagate_inbounds function _sqrt_quasitriu_offdiag_block_1x1!(R, A, i, j)
2667
    Rii = R[i, i]
65,734✔
2668
    Rjj = R[j, j]
65,734✔
2669
    iszero(Rii) && iszero(Rjj) && return R
65,734✔
2670
    t = eltype(R)
65,733✔
2671
    tt = typeof(zero(t)*zero(t))
65,733✔
2672
    r = tt(-A[i, j])
65,733✔
2673
    @simd for k in (i + 1):(j - 1)
72,656✔
2674
        r += R[i, k] * R[k, j]
2,991,886✔
2675
    end
2676
    iszero(r) && return R
65,733✔
2677
    R[i, j] = sylvester(Rii, Rjj, r)
56,712✔
2678
    return R
56,712✔
2679
end
2680

2681
Base.@propagate_inbounds function _sqrt_quasitriu_offdiag_block_1x2!(R, A, i, j)
2682
    jrange = j:(j + 1)
250✔
2683
    t = eltype(R)
250✔
2684
    tt = typeof(zero(t)*zero(t))
250✔
2685
    r1 = tt(-A[i, j])
250✔
2686
    r2 = tt(-A[i, j + 1])
250✔
2687
    @simd for k in (i + 1):(j - 1)
360✔
2688
        rik = R[i, k]
476✔
2689
        r1 += rik * R[k, j]
476✔
2690
        r2 += rik * R[k, j + 1]
476✔
2691
    end
2692
    Rjj = @view R[jrange, jrange]
250✔
2693
    Rij = @view R[i, jrange]
250✔
2694
    Rij[1] = r1
250✔
2695
    Rij[2] = r2
250✔
2696
    _sylvester_1x2!(R[i, i], Rjj, Rij)
250✔
2697
    return R
250✔
2698
end
2699

2700
Base.@propagate_inbounds function _sqrt_quasitriu_offdiag_block_2x1!(R, A, i, j)
2701
    irange = i:(i + 1)
417✔
2702
    t = eltype(R)
417✔
2703
    tt = typeof(zero(t)*zero(t))
417✔
2704
    r1 = tt(-A[i, j])
417✔
2705
    r2 = tt(-A[i + 1, j])
417✔
2706
    @simd for k in (i + 2):(j - 1)
546✔
2707
        rkj = R[k, j]
1,003✔
2708
        r1 += R[i, k] * rkj
1,003✔
2709
        r2 += R[i + 1, k] * rkj
1,003✔
2710
    end
2711
    Rii = @view R[irange, irange]
417✔
2712
    Rij = @view R[irange, j]
417✔
2713
    Rij[1] = r1
417✔
2714
    Rij[2] = r2
417✔
2715
    @views _sylvester_2x1!(Rii, R[j, j], Rij)
417✔
2716
    return R
417✔
2717
end
2718

2719
Base.@propagate_inbounds function _sqrt_quasitriu_offdiag_block_2x2!(R, A, i, j)
2720
    irange = i:(i + 1)
414✔
2721
    jrange = j:(j + 1)
414✔
2722
    t = eltype(R)
414✔
2723
    tt = typeof(zero(t)*zero(t))
414✔
2724
    for i′ in irange, j′ in jrange
414✔
2725
        Cij = tt(-A[i′, j′])
1,656✔
2726
        @simd for k in (i + 2):(j - 1)
2,332✔
2727
            Cij += R[i′, k] * R[k, j′]
3,464✔
2728
        end
2729
        R[i′, j′] = Cij
1,656✔
2730
    end
2,070✔
2731
    Rii = @view R[irange, irange]
414✔
2732
    Rjj = @view R[jrange, jrange]
414✔
2733
    Rij = @view R[irange, jrange]
414✔
2734
    if !iszero(Rij) && !all(isnan, Rij)
417✔
2735
        _sylvester_2x2!(Rii, Rjj, Rij)
411✔
2736
    end
2737
    return R
414✔
2738
end
2739

2740
# solve Sylvester's equation AX + XB = -C using blockwise recursion until the dimension of
2741
# A and B are no greater than blockwidth, based on Algorithm 1 from
2742
# Jonsson I, Kågström B. Recursive blocked algorithms for solving triangular systems—
2743
# Part I: one-sided and coupled Sylvester-type matrix equations. (2002) ACM Trans Math Softw.
2744
# 28(4), https://doi.org/10.1145/592843.592845.
2745
# specify raise=false to avoid breaking the recursion if a LAPACKException is thrown when
2746
# computing one of the blocks.
2747
function _sylvester_quasitriu!(A, B, C; blockwidth=64, nA=checksquare(A), nB=checksquare(B), raise=true)
863✔
2748
    if 1 ≤ nA ≤ blockwidth && 1 ≤ nB ≤ blockwidth
539✔
2749
        _sylvester_quasitriu_base!(A, B, C; raise=raise)
398✔
2750
    elseif nA ≥ 2nB ≥ 2
141✔
2751
        _sylvester_quasitriu_split1!(A, B, C; blockwidth=blockwidth, nA=nA, nB=nB, raise=raise)
12✔
2752
    elseif nB ≥ 2nA ≥ 2
129✔
2753
        _sylvester_quasitriu_split2!(A, B, C; blockwidth=blockwidth, nA=nA, nB=nB, raise=raise)
12✔
2754
    else
2755
        _sylvester_quasitriu_splitall!(A, B, C; blockwidth=blockwidth, nA=nA, nB=nB, raise=raise)
117✔
2756
    end
2757
    return C
499✔
2758
end
2759
function _sylvester_quasitriu_base!(A, B, C; raise=true)
796✔
2760
    try
398✔
2761
        _, scale = LAPACK.trsyl!('N', 'N', A, B, C)
398✔
2762
        rmul!(C, -inv(scale))
398✔
2763
    catch e
2764
        if !(e isa LAPACKException) || raise
296✔
2765
            throw(e)
148✔
2766
        end
2767
    end
2768
    return C
382✔
2769
end
2770
function _sylvester_quasitriu_split1!(A, B, C; nA=checksquare(A), kwargs...)
24✔
2771
    iA = div(nA, 2)
12✔
2772
    iszero(A[iA + 1, iA]) || (iA += 1)  # don't split 2x2 diagonal block
12✔
2773
    rA1, rA2 = 1:iA, (iA + 1):nA
12✔
2774
    nA1, nA2 = iA, nA-iA
12✔
2775
    A11, A12, A22 = @views A[rA1,rA1], A[rA1,rA2], A[rA2,rA2]
12✔
2776
    C1, C2 = @views C[rA1,:], C[rA2,:]
12✔
2777
    _sylvester_quasitriu!(A22, B, C2; nA=nA2, kwargs...)
24✔
2778
    mul!(C1, A12, C2, true, true)
8✔
2779
    _sylvester_quasitriu!(A11, B, C1; nA=nA1, kwargs...)
16✔
2780
    return C
8✔
2781
end
2782
function _sylvester_quasitriu_split2!(A, B, C; nB=checksquare(B), kwargs...)
24✔
2783
    iB = div(nB, 2)
12✔
2784
    iszero(B[iB + 1, iB]) || (iB += 1)  # don't split 2x2 diagonal block
12✔
2785
    rB1, rB2 = 1:iB, (iB + 1):nB
12✔
2786
    nB1, nB2 = iB, nB-iB
12✔
2787
    B11, B12, B22 = @views B[rB1,rB1], B[rB1,rB2], B[rB2,rB2]
12✔
2788
    C1, C2 = @views C[:,rB1], C[:,rB2]
12✔
2789
    _sylvester_quasitriu!(A, B11, C1; nB=nB1, kwargs...)
24✔
2790
    mul!(C2, C1, B12, true, true)
8✔
2791
    _sylvester_quasitriu!(A, B22, C2; nB=nB2, kwargs...)
16✔
2792
    return C
8✔
2793
end
2794
function _sylvester_quasitriu_splitall!(A, B, C; nA=checksquare(A), nB=checksquare(B), kwargs...)
234✔
2795
    iA = div(nA, 2)
117✔
2796
    iszero(A[iA + 1, iA]) || (iA += 1)  # don't split 2x2 diagonal block
128✔
2797
    iB = div(nB, 2)
117✔
2798
    iszero(B[iB + 1, iB]) || (iB += 1)  # don't split 2x2 diagonal block
136✔
2799
    rA1, rA2 = 1:iA, (iA + 1):nA
117✔
2800
    nA1, nA2 = iA, nA-iA
117✔
2801
    rB1, rB2 = 1:iB, (iB + 1):nB
117✔
2802
    nB1, nB2 = iB, nB-iB
117✔
2803
    A11, A12, A22 = @views A[rA1,rA1], A[rA1,rA2], A[rA2,rA2]
117✔
2804
    B11, B12, B22 = @views B[rB1,rB1], B[rB1,rB2], B[rB2,rB2]
117✔
2805
    C11, C21, C12, C22 = @views C[rA1,rB1], C[rA2,rB1], C[rA1,rB2], C[rA2,rB2]
117✔
2806
    _sylvester_quasitriu!(A22, B11, C21; nA=nA2, nB=nB1, kwargs...)
129✔
2807
    mul!(C11, A12, C21, true, true)
101✔
2808
    _sylvester_quasitriu!(A11, B11, C11; nA=nA1, nB=nB1, kwargs...)
109✔
2809
    mul!(C22, C21, B12, true, true)
101✔
2810
    _sylvester_quasitriu!(A22, B22, C22; nA=nA2, nB=nB2, kwargs...)
109✔
2811
    mul!(C12, A12, C22, true, true)
101✔
2812
    mul!(C12, C11, B12, true, true)
101✔
2813
    _sylvester_quasitriu!(A11, B22, C12; nA=nA1, nB=nB2, kwargs...)
109✔
2814
    return C
101✔
2815
end
2816

2817
# End of auxiliary functions for matrix square root
2818

2819
# Generic eigensystems
2820
eigvals(A::AbstractTriangular) = diag(A)
108✔
2821
# fallback for unknown types
2822
function eigvecs(A::AbstractTriangular{<:BlasFloat})
6✔
2823
    if istriu(A)
8✔
2824
        eigvecs(UpperTriangular(Matrix(A)))
4✔
2825
    else # istril(A)
2826
        eigvecs(LowerTriangular(Matrix(A)))
4✔
2827
    end
2828
end
2829
function eigvecs(A::AbstractTriangular{T}) where T
8✔
2830
    TT = promote_type(T, Float32)
12✔
2831
    if TT <: BlasFloat
12✔
2832
        return eigvecs(convert(AbstractMatrix{TT}, A))
12✔
2833
    else
2834
        throw(ArgumentError(lazy"eigvecs type $(typeof(A)) not supported. Please submit a pull request."))
×
2835
    end
2836
end
2837
det(A::UnitUpperTriangular{T}) where {T} = one(T)
7✔
2838
det(A::UnitLowerTriangular{T}) where {T} = one(T)
7✔
2839
logdet(A::UnitUpperTriangular{T}) where {T} = zero(T)
7✔
2840
logdet(A::UnitLowerTriangular{T}) where {T} = zero(T)
7✔
2841
logabsdet(A::UnitUpperTriangular{T}) where {T} = zero(T), one(T)
7✔
2842
logabsdet(A::UnitLowerTriangular{T}) where {T} = zero(T), one(T)
7✔
2843
det(A::UpperTriangular) = prod(diag(A.data))
308✔
2844
det(A::LowerTriangular) = prod(diag(A.data))
7✔
2845
function logabsdet(A::Union{UpperTriangular{T},LowerTriangular{T}}) where T
65✔
2846
    sgn = one(T)
104✔
2847
    abs_det = zero(real(T))
104✔
2848
    @inbounds for i in 1:size(A,1)
104✔
2849
        diag_i = A.data[i,i]
664✔
2850
        sgn *= sign(diag_i)
826✔
2851
        abs_det += log(abs(diag_i))
775✔
2852
    end
1,224✔
2853
    return abs_det, sgn
104✔
2854
end
2855

2856
eigen(A::AbstractTriangular) = Eigen(eigvals(A), eigvecs(A))
20✔
2857

2858
# Generic singular systems
2859
for func in (:svd, :svd!, :svdvals)
2860
    @eval begin
2861
        ($func)(A::AbstractTriangular; kwargs...) = ($func)(copyto!(similar(parent(A)), A); kwargs...)
112✔
2862
    end
2863
end
2864

2865
factorize(A::AbstractTriangular) = A
28✔
2866

2867
# disambiguation methods: /(Adjoint of AbsVec, <:AbstractTriangular)
2868
/(u::AdjointAbsVec, A::Union{LowerTriangular,UpperTriangular}) = adjoint(adjoint(A) \ u.parent)
×
2869
/(u::AdjointAbsVec, A::Union{UnitLowerTriangular,UnitUpperTriangular}) = adjoint(adjoint(A) \ u.parent)
×
2870
# disambiguation methods: /(Transpose of AbsVec, <:AbstractTriangular)
2871
/(u::TransposeAbsVec, A::Union{LowerTriangular,UpperTriangular}) = transpose(transpose(A) \ u.parent)
×
2872
/(u::TransposeAbsVec, A::Union{UnitLowerTriangular,UnitUpperTriangular}) = transpose(transpose(A) \ u.parent)
×
2873
# disambiguation methods: /(Transpose of AbsVec, Adj/Trans of <:AbstractTriangular)
2874
for (tritype, comptritype) in ((:LowerTriangular, :UpperTriangular),
2875
                               (:UnitLowerTriangular, :UnitUpperTriangular),
2876
                               (:UpperTriangular, :LowerTriangular),
2877
                               (:UnitUpperTriangular, :UnitLowerTriangular))
2878
    @eval /(u::TransposeAbsVec, A::$tritype{<:Any,<:Adjoint}) = transpose($comptritype(conj(parent(parent(A)))) \ u.parent)
×
2879
    @eval /(u::TransposeAbsVec, A::$tritype{<:Any,<:Transpose}) = transpose(transpose(A) \ u.parent)
×
2880
end
2881

2882
# Cube root of a 2x2 real-valued matrix with complex conjugate eigenvalues and equal diagonal values.
2883
# Reference [1]: Smith, M. I. (2003). A Schur Algorithm for Computing Matrix pth Roots.
2884
#   SIAM Journal on Matrix Analysis and Applications (Vol. 24, Issue 4, pp. 971–989).
2885
#   https://doi.org/10.1137/s0895479801392697
2886
function _cbrt_2x2!(A::AbstractMatrix{T}) where {T<:Real}
3✔
2887
    @assert checksquare(A) == 2
3✔
2888
    @inbounds begin
3✔
2889
        (A[1,1] == A[2,2]) || throw(ArgumentError("_cbrt_2x2!: Matrix A must have equal diagonal values."))
3✔
2890
        (A[1,2]*A[2,1] < 0) || throw(ArgumentError("_cbrt_2x2!: Matrix A must have complex conjugate eigenvalues."))
3✔
2891
        μ = sqrt(-A[1,2]*A[2,1])
3✔
2892
        r = cbrt(hypot(A[1,1], μ))
3✔
2893
        θ = atan(μ, A[1,1])
3✔
2894
        s, c = sincos(θ/3)
3✔
2895
        α, β′ = r*c, r*s/µ
3✔
2896
        A[1,1] = α
3✔
2897
        A[2,2] = α
3✔
2898
        A[1,2] = β′*A[1,2]
3✔
2899
        A[2,1] = β′*A[2,1]
3✔
2900
    end
2901
    return A
3✔
2902
end
2903

2904
# Cube root of a quasi upper triangular matrix (output of Schur decomposition)
2905
# Reference [1]: Smith, M. I. (2003). A Schur Algorithm for Computing Matrix pth Roots.
2906
#   SIAM Journal on Matrix Analysis and Applications (Vol. 24, Issue 4, pp. 971–989).
2907
#   https://doi.org/10.1137/s0895479801392697
2908
@views function _cbrt_quasi_triu!(A::AbstractMatrix{T}) where {T<:Real}
3✔
2909
    m, n = size(A)
3✔
2910
    (m == n) || throw(ArgumentError("_cbrt_quasi_triu!: Matrix A must be square."))
3✔
2911
    # Cube roots of 1x1 and 2x2 diagonal blocks
2912
    i = 1
3✔
2913
    sizes = ones(Int,n)
30✔
2914
    S = zeros(T,2,n)
60✔
2915
    while i < n
27✔
2916
        if !iszero(A[i+1,i])
24✔
2917
            _cbrt_2x2!(A[i:i+1,i:i+1])
3✔
2918
            mul!(S[1:2,i:i+1], A[i:i+1,i:i+1], A[i:i+1,i:i+1])
3✔
2919
            sizes[i] = 2
3✔
2920
            sizes[i+1] = 0
3✔
2921
            i += 2
3✔
2922
        else
2923
            A[i,i] = cbrt(A[i,i])
21✔
2924
            S[1,i] = A[i,i]*A[i,i]
21✔
2925
            i += 1
21✔
2926
        end
2927
    end
24✔
2928
    if i == n
3✔
2929
        A[n,n] = cbrt(A[n,n])
3✔
2930
        S[1,n] = A[n,n]*A[n,n]
3✔
2931
    end
2932
    # Algorithm 4.3 in Reference [1]
2933
    Δ = I(4)
12✔
2934
    M_L₀ = zeros(T,4,4)
48✔
2935
    M_L₁ = zeros(T,4,4)
48✔
2936
    M_Bᵢⱼ⁽⁰⁾ = zeros(T,2,2)
12✔
2937
    M_Bᵢⱼ⁽¹⁾ = zeros(T,2,2)
12✔
2938
    for k = 1:n-1
3✔
2939
        for i = 1:n-k
27✔
2940
            if sizes[i] == 0 || sizes[i+k] == 0 continue end
255✔
2941
            k₁, k₂ = i+1+(sizes[i+1]==0), i+k-1
111✔
2942
            i₁, i₂, j₁, j₂, s₁, s₂ = i, i+sizes[i]-1, i+k, i+k+sizes[i+k]-1, sizes[i], sizes[i+k]
111✔
2943
            L₀ = M_L₀[1:s₁*s₂,1:s₁*s₂]
111✔
2944
            L₁ = M_L₁[1:s₁*s₂,1:s₁*s₂]
111✔
2945
            Bᵢⱼ⁽⁰⁾ = M_Bᵢⱼ⁽⁰⁾[1:s₁, 1:s₂]
111✔
2946
            Bᵢⱼ⁽¹⁾ = M_Bᵢⱼ⁽¹⁾[1:s₁, 1:s₂]
111✔
2947
            # Compute Bᵢⱼ⁽⁰⁾ and Bᵢⱼ⁽¹⁾
2948
            mul!(Bᵢⱼ⁽⁰⁾, A[i₁:i₂,k₁:k₂], A[k₁:k₂,j₁:j₂])
111✔
2949
            # Retrieve Rᵢ,ᵢ₊ₖ as A[i+k,i]'
2950
            mul!(Bᵢⱼ⁽¹⁾, A[i₁:i₂,k₁:k₂], A[j₁:j₂,k₁:k₂]')
111✔
2951
            # Solve Uᵢ,ᵢ₊ₖ using Reference [1, (4.10)]
2952
            kron!(L₀, Δ[1:s₂,1:s₂], S[1:s₁,i₁:i₂])
111✔
2953
            L₀ .+= kron!(L₁, A[j₁:j₂,j₁:j₂]', A[i₁:i₂,i₁:i₂])
222✔
2954
            L₀ .+= kron!(L₁, S[1:s₂,j₁:j₂]', Δ[1:s₁,1:s₁])
222✔
2955
            mul!(A[i₁:i₂,j₁:j₂], A[i₁:i₂,i₁:i₂], Bᵢⱼ⁽⁰⁾, -1.0, 1.0)
111✔
2956
            A[i₁:i₂,j₁:j₂] .-= Bᵢⱼ⁽¹⁾
222✔
2957
            ldiv!(lu!(L₀), A[i₁:i₂,j₁:j₂][:])
111✔
2958
            # Compute and store Rᵢ,ᵢ₊ₖ' in A[i+k,i]
2959
            mul!(Bᵢⱼ⁽⁰⁾, A[i₁:i₂,i₁:i₂], A[i₁:i₂,j₁:j₂], 1.0, 1.0)
111✔
2960
            mul!(Bᵢⱼ⁽⁰⁾, A[i₁:i₂,j₁:j₂], A[j₁:j₂,j₁:j₂], 1.0, 1.0)
111✔
2961
            A[j₁:j₂,i₁:i₂] .= Bᵢⱼ⁽⁰⁾'
111✔
2962
        end
243✔
2963
    end
51✔
2964
    # Make quasi triangular
2965
    for j=1:n for i=j+1+(sizes[j]==2):n A[i,j] = 0 end end
30✔
2966
    return A
3✔
2967
end
2968

2969
# Cube roots of real-valued triangular matrices
2970
cbrt(A::UpperTriangular{T}) where {T<:Real} = UpperTriangular(_cbrt_quasi_triu!(Matrix{T}(A)))
1✔
2971
cbrt(A::LowerTriangular{T}) where {T<:Real} = LowerTriangular(_cbrt_quasi_triu!(Matrix{T}(A'))')
1✔
STATUS · Troubleshooting · Open an Issue · Sales · Support · CAREERS · ENTERPRISE · START FREE · SCHEDULE DEMO
ANNOUNCEMENTS · TWITTER · TOS & SLA · Supported CI Services · What's a CI service? · Automated Testing

© 2026 Coveralls, Inc