-
Notifications
You must be signed in to change notification settings - Fork 3
/
tvdiff.jl
316 lines (277 loc) · 9.79 KB
/
tvdiff.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
"""
tvdiff(data::AbstractVector, iter::Integer, α::Real; kwargs...)
# Arguments
- `data::AbstractVector`:
Vector of data to be differentiated.
- `iter::Integer`:
Number of iterations to run the main loop. A stopping
condition based on the norm of the gradient vector `g`
below would be an easy modification.
- `α::Real`:
Regularization parameter. This is the main parameter
to fiddle with. Start by varying by orders of
magnitude until reasonable results are obtained. A
value to the nearest power of 10 is usally adequate.
Higher values increase regularization strenght
and improve conditioning.
## Keywords
- `u_0::AbstractVector`:
Initialization of the iteration. Default value is the
naive derivative (without scaling), of appropriate
length (this being different for the two methods).
Although the solution is theoretically independent of
the intialization, a poor choice can exacerbate
conditioning issues when the linear system is solved.
- `scale::String`:
Scale of dataset, `\"large\"` or `\"small\"` (case insensitive).
Default is `\"small\"` . `\"small\"` has somewhat better
boundary behavior, but becomes unwieldly for very large datasets.
`\"large\"` has simpler numerics but
is more efficient for large-scale problems. `\"large\"` is
more readily modified for higher-order derivatives,
since the implicit differentiation matrix is square.
- `ε::Real`:
Parameter for avoiding division by zero. Default value
is `1e-6`. Results should not be very sensitive to the
value. Larger values improve conditioning and
therefore speed, while smaller values give more
accurate results with sharper jumps.
- `dx::Real`:
Grid spacing, used in the definition of the derivative
operators. Default is `1 / length(data)`.
- `precond::String`:
Select the preconditioner for the conjugate gradient method.
Default is `\"none\"`.
+ `scale = \"small\"`:
While in principle `precond=\"simple\"` should speed things up,
sometimes the preconditioner can cause convergence problems instead,
and should be left to `\"none\"`.
+ `scale = \"large\"`:
The improved preconditioners are one of the main features of the
algorithm, therefore using the default `\"none\"` is discouraged.
Currently, `\"diagonal\"`,`\"amg_rs\"`,`\"amg_sa\"`, `\"cholesky\"` are available.
- `diff_kernel::String`:
Kernel to use in the integral to smooth the derivative. By default it is set to
`\"abs\"`, the absolute value ``|u'|``. However, it can be changed to `\"square\"`,
the square value ``(u')^2``. The latter produces smoother
derivatives, whereas the absolute values tends to make them more blocky.
- `cg_tol::Real`:
Relative tolerance used in conjugate gradient method. Default is `1e-6`.
- `cg_maxiter::Int`:
Maximum number of iterations to use in conjugate gradient optimisation.
Default is `100`.
- `show_diagn::Bool`:
Flag whether to display diagnostics at each iteration. Default is `false`.
Useful for diagnosing preconditioning problems. When tolerance is not met,
an early iterate being best is more worrying than a large relative residual.
# Output
- `u`:
Estimate of the regularized derivative of data with
`length(u) = length(data)`.
"""
function tvdiff(
data::AbstractVector,
iter::Integer,
α::Real;
u_0::AbstractVector=[NaN],
scale::String="small",
ε::Real=1e-6,
dx::Real=NaN,
precond::String="none",
diff_kernel::String="abs",
cg_tol::Real=1e-6,
cg_maxiter::Integer=100,
show_diagn::Bool=false,
)::AbstractVector
n = length(data)
if isnan(dx)
dx = 1 / (n - 1)
end
# Make string inputs case insensitive
scale = lowercase(scale)
precond = lowercase(precond)
diff_kernel = lowercase(diff_kernel)
# Run tvdiff for selected method
return tvdiff(
Val(Symbol(scale)),
data,
iter,
α,
u_0,
ε,
dx,
cg_tol,
cg_maxiter,
precond,
diff_kernel,
show_diagn,
)
end
# Total variance regularized numerical differences for small scale problems.
function tvdiff(
::Val{:small},
data::AbstractVector,
iter::Integer,
α::Real,
u_0::AbstractVector,
ε::Real,
dx::Real,
cg_tol::Real,
cg_maxiter::Integer,
precond::String,
diff_kernel::String,
show_diagn::Bool,
)::AbstractVector
n = length(data)
#= Assert initialization if provided, otherwise set
default initization to naive derivative =#
if isequal(u_0, [NaN])
u_0 = [0; diff(data); 0] / dx
elseif length(u_0) != (n + 1)
throw(
DimensionMismatch(
"size $(size(u_0)) of u_0 doesn't match size ($(n + 1),) required for scale=\"small\".",
),
)
end
u = copy(u_0)
# Construct differentiation matrix.
D = spdiagm(n, n + 1, 0 => -ones(n), 1 => ones(n)) / dx
Dᵀ = transpose(D)
# Construct antidifferentiation operator and its adjoint.
function A(x)
return (cumsum(x) - 0.5 * (x .- x[1]))[2:end] * dx
end
function Aᵀ(x)
return [sum(x) / 2; (sum(x) .- cumsum(x) .- x / 2)] * dx
end
# Precompute antidifferentiation adjoint on data
# Since A([0]) = 0, we need to adjust.
offset = data[1]
Aᵀb = Aᵀ(offset .- data)
for i in 1:iter
if diff_kernel == "abs"
# Diagonal matrix of weights, for linearizing E-L equation.
Q = Diagonal(1 ./ sqrt.((D * u) .^ 2 .+ ε))
# Linearized diffusion matrix, also approximation of Hessian.
L = dx * Dᵀ * Q * D
elseif diff_kernel == "square"
L = dx * Dᵀ * D
else
throw(
ArgumentError(
"""unexpected diff_kernel "$(diff_kernel)" for scale="small"."""
),
)
end
# Gradient of functional.
g = Aᵀ(A(u)) + Aᵀb + α * L * u
# Select preconditioner.
if precond == "simple"
P = Diagonal(α * diag(L) .+ 1)
elseif precond == "none"
P = I # Identity matrix
else
throw(ArgumentError("""unexpected precond "$(precond)" for scale="small"."""))
end
# Prepare linear operator for linear equation.
# Approximation of Hessian of TVR functional at u
H = LinearMap(u -> Aᵀ(A(u)) + α * L * u, n + 1, n + 1)
# Solve linear equation.
s = cg(H, -g; Pl=P, reltol=cg_tol, maxiter=cg_maxiter)
show_diagn && println(
"Iteration $(i):\trel. change = $(norm(s) / norm(u)),\tgradient norm = $(norm(g))",
)
# Update current solution
u += s
end
return u[1:(end - 1)]
end
# Total variance regularized numerical differences for large scale problems.
function tvdiff(
::Val{:large},
data::AbstractVector,
iter::Integer,
α::Real,
u_0::AbstractVector,
ε::Real,
dx::Real,
cg_tol::Real,
cg_maxiter::Integer,
precond::String,
diff_kernel::String,
show_diagn::Bool,
)::AbstractVector
n = length(data)
# Since Au(0) = 0, we need to adjust.
data = data .- data[1]
#= Assert initialization if provided, otherwise set
default initization to naive derivative =#
if isequal(u_0, [NaN])
u_0 = [0; diff(data)] / dx
elseif length(u_0) != n
throw(
DimensionMismatch(
"size $(size(u_0)) of u_0 doesn't match size ($(n),) required for scale=\"large\".",
),
)
end
u = deepcopy(u_0) * dx
# Construct differentiation matrix.
D = spdiagm(n, n, 0 => -ones(n - 1), 1 => ones(n - 1)) / dx
Dᵀ = transpose(D)
# Construct antidifferentiation operator and its adjoint.
A(x) = cumsum(x)
function Aᵀ(x)
return sum(x) .- [0; cumsum(x[1:(end - 1)])]
end
# Precompute antidifferentiation adjoint on data
Aᵀd = Aᵀ(data)
for i in 1:iter
if diff_kernel == "abs"
# Diagonal matrix of weights, for linearizing E-L equation.
Q = Diagonal(1 ./ sqrt.((D * u) .^ 2 .+ ε))
# Linearized diffusion matrix, also approximation of Hessian.
L = Dᵀ * Q * D
elseif diff_kernel == "square"
L = Dᵀ * D
else
throw(
ArgumentError(
"""unexpected diff_kernel "$(diff_kernel)" for scale="large"."""
),
)
end
# Gradient of functional.
g = Aᵀ(A(u)) - Aᵀd + α * L * u
# Select preconditioner.
B = α * L + Diagonal(reverse(cumsum(n:-1:1)))
if precond == "cholesky"
# Incomplete Cholesky preconditioner with cut-off level 2
P = CholeskyPreconditioner(B, 2)
elseif precond == "diagonal"
P = DiagonalPreconditioner(B)
elseif precond == "amg_rs"
# Ruge-Stuben variant
P = AMGPreconditioner{RugeStuben}(B)
elseif precond == "amg_sa"
# Smoothed aggregation
P = AMGPreconditioner{SmoothedAggregation}(B)
elseif precond == "none"
P = I # Identity matrix
else
throw(ArgumentError("""unexpected precod "$(precond)" for scale="large"."""))
end
# Prepare linear operator for linear equation.
# Approximation of Hessian of TVR functional at u
H = LinearMap(u -> Aᵀ(A(u)) + α * L * u, n, n)
# Solve linear equation.
s = cg(H, -g; Pl=P, reltol=cg_tol, maxiter=cg_maxiter)
show_diagn && println(
"Iteration $(i):\trel. change = $(norm(s) / norm(u)),\tgradient norm = $(norm(g))",
)
# Update current solution
u += s
end
return u / dx
end