-
Notifications
You must be signed in to change notification settings - Fork 1
/
DifferentiableFactorizations.jl
198 lines (178 loc) · 4.77 KB
/
DifferentiableFactorizations.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
module DifferentiableFactorizations
export diff_qr, diff_cholesky, diff_lu, diff_eigen, diff_svd, diff_schur
using LinearAlgebra, ImplicitDifferentiation, ComponentArrays, ChainRulesCore
# QR
function qr_conditions(A, x)
(; Q, R) = x
return vcat(
vec(UpperTriangular(Q' * Q) + LowerTriangular(R) - I - Diagonal(R)),
vec(Q * R - A),
)
end
function qr_forward(A)
qr_res = qr(A)
Q = copy(qr_res.Q[:, 1:size(A, 2)])
(; R) = qr_res
return ComponentVector(; Q, R)
end
const _diff_qr = ImplicitFunction(qr_forward, qr_conditions, DirectLinearSolver(), nothing)
function diff_qr(A)
(; Q, R) = _diff_qr(A)
return (; Q, R)
end
# Cholesky
function cholesky_conditions(A, U)
return vec(
UpperTriangular(U' * U) + LowerTriangular(U) - UpperTriangular(A) - Diagonal(U),
)
end
function cholesky_forward(A)
ch_res = cholesky(A)
return ch_res.U
end
const _diff_cholesky =
ImplicitFunction(cholesky_forward, cholesky_conditions, DirectLinearSolver(), nothing)
function diff_cholesky(A)
U = _diff_cholesky(A)
return (; L = U', U)
end
# LU
function lu_conditions(A, LU, p)
(; L, U) = LU
return vcat(
vec(UpperTriangular(L) + LowerTriangular(U) - Diagonal(U) - I),
vec(L * U - A[p, :]),
)
end
function lu_forward(A)
lu_res = lu(A)
(; L, U, p) = lu_res
return ComponentVector(; L, U), p
end
const _diff_lu = ImplicitFunction(lu_forward, lu_conditions, DirectLinearSolver(), nothing)
function diff_lu(A)
temp, p = _diff_lu(A)
(; L, U) = temp
return (; L, U, p)
end
# Eigen
comp_vec(A) = ComponentVector((; A))
comp_vec(A, B) = ComponentVector((; A, B))
function ChainRulesCore.rrule(::typeof(comp_vec), A)
out = comp_vec(A)
T = typeof(out)
return out, Δ -> begin
_Δ = convert(T, Δ)
(NoTangent(), _Δ.A)
end
end
function ChainRulesCore.rrule(::typeof(comp_vec), A, B)
out = comp_vec(A, B)
T = typeof(out)
return out, Δ -> begin
_Δ = convert(T, Δ)
(NoTangent(), _Δ.A, _Δ.B)
end
end
function eigen_conditions(AB, sV)
(; s, V) = sV
(; A) = AB
if hasproperty(AB, :B)
(; B) = AB
else
B = I
end
return vcat(vec(A * V - B * V * Diagonal(s)), diag(V' * B * V) .- 1)
end
function eigen_forward(AB)
(; A) = AB
if hasproperty(AB, :B)
(; B) = AB
eig_res = eigen(A, B)
else
eig_res = eigen(A)
end
s = eig_res.values
V = eig_res.vectors
return ComponentVector(; s, V)
end
const _diff_eigen =
ImplicitFunction(eigen_forward, eigen_conditions, DirectLinearSolver(), nothing)
function diff_eigen(A)
(; s, V) = _diff_eigen(comp_vec(A))
return (; s, V)
end
function diff_eigen(A, B)
(; s, V) = _diff_eigen(comp_vec(A, B))
return (; s, V)
end
function schur_conditions(A, Z_T)
(; Z, T) = Z_T
return vcat(vec(Z' * A * Z - T), vec(Z' * Z - I + LowerTriangular(T) - Diagonal(T)))
end
function schur_forward(A)
schur_res = schur(A)
(; Z, T) = schur_res
return ComponentVector(; Z, T)
end
const _diff_schur =
ImplicitFunction(schur_forward, schur_conditions, DirectLinearSolver(), nothing)
function bidiag(v1, v2)
return Bidiagonal(v1, v2, :L)
end
function ChainRulesCore.rrule(::typeof(bidiag), v1, v2)
bidiag(v1, v2), Δ -> begin
NoTangent(), diag(Δ), diag(Δ, -1)
end
end
function gen_schur_conditions(AB, left_right_S_T)
(; left, right, S, T) = left_right_S_T
(; A, B) = AB
return vcat(
vec(left * S * right' - A),
vec(left * T * right' - B),
vec(
UpperTriangular(left' * left) - I + LowerTriangular(S) -
bidiag(diag(S), diag(S, -1) .+ (diag(S, -1) .* diag(T, 1))),
),
vec(UpperTriangular(right' * right) - I + LowerTriangular(T) - Diagonal(T)),
)
end
function gen_schur_forward(AB)
(; A, B) = AB
schur_res = schur(A, B)
(; left, right, S, T) = schur_res
return ComponentVector(; left, right, S, T)
end
const _diff_gen_schur =
ImplicitFunction(gen_schur_forward, gen_schur_conditions, DirectLinearSolver(), nothing)
function diff_schur(A, B)
(; left, right, S, T) = _diff_gen_schur(comp_vec(A, B))
return (; left, right, S, T)
end
function diff_schur(A)
(; Z, T) = _diff_schur(A)
return (; Z, T)
end
# SVD
function svd_conditions(A, USV)
(; U, S, V) = USV
VtV = V' * V
return vcat(
vec(U * Diagonal(S) * V' - A),
vec(UpperTriangular(VtV) + LowerTriangular(U' * U) - 2I),
diag(VtV) .- 1,
)
end
function svd_forward(A)
svd_res = svd(A)
(; U, S, V) = svd_res
return ComponentVector(; U, S, V)
end
const _diff_svd =
ImplicitFunction(svd_forward, svd_conditions, DirectLinearSolver(), nothing)
function diff_svd(A)
(; U, S, V) = _diff_svd(A)
return (; U, S, V)
end
end