forked from jeremylt/LFAToolkit.jl
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathUtilities.jl
372 lines (313 loc) · 10.6 KB
/
Utilities.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
# ------------------------------------------------------------------------------
# convenience utilities
# ------------------------------------------------------------------------------
# ------------------------------------------------------------------------------
# compute symbol matrix
# ------------------------------------------------------------------------------
"""
```julia
computesymbolsoverrange(operator, θ; mass=nothing, θ_min=-π/2)
```
Compute the eigenvalues and eigenvectors of the symbol matrix for an operator over
a range of θ from -π/2 to 3π/2
# Arguments:
- `operator`: Finite element operator to compute symbol matrices for
- `numbersteps1d`: Number of values of θ to sample in each dimension
Note: numbersteps1d^dimension symbol matrices will be computed
# Keyword Arguments:
- `mass`: Mass operator to invert for comparison to analytic solution,
default: nothing
- `θ_min`: Bottom of range of θ, shifts range to [θ_min, θ_min + 2π],
default: -π / 2
# Returns:
- Values of θ sampled
- Eigenvalues of symbol matrix at θ sampled
- Eigenvectors of symbol matrix at θ sampled
# Examples:
```jldoctest
using LFAToolkit;
numbersteps1d = 5
for dimension in 1:3
# setup
mesh = []
if dimension == 1
mesh = Mesh1D(1.0);
elseif dimension == 2
mesh = Mesh2D(1.0, 1.0);
elseif dimension == 3
mesh = Mesh3D(1.0, 1.0, 1.0);
end
diffusion = GalleryOperator("diffusion", 3, 3, mesh);
# compute symbols
(_, eigenvalues, _) = computesymbolsoverrange(diffusion, numbersteps1d);
# verify
eigenvalues = real(eigenvalues);
if dimension == 1
@assert minimum(eigenvalues[4, :]) ≈ 1
@assert maximum(eigenvalues[4, :]) ≈ 4/3
elseif dimension == 2
@assert minimum(eigenvalues[19, :]) ≈ 2/3
@assert maximum(eigenvalues[19, :]) ≈ 64/45
elseif dimension == 3
@assert minimum(eigenvalues[94, :]) ≈ 1/3
@assert maximum(eigenvalues[94, :]) ≈ 256/225
end
end
# output
```
```jldoctest
using LFAToolkit;
# setup
numbersteps1d = 3
mesh = Mesh2D(0.5, 0.5)
p = 6
# operators
mass = GalleryOperator("mass", p + 1, p + 1, mesh)
diffusion = GalleryOperator("diffusion", p + 1, p + 1, mesh)
# compute symbols
(_, eigenvalues, _) = computesymbolsoverrange(
diffusion,
numbersteps1d;
mass = mass,
θ_min = -π
);
eigenvalues = real(eigenvalues);
@assert minimum(eigenvalues[2, :]) ≈ π^2
# output
```
"""
function computesymbolsoverrange(
operator::Operator,
numbersteps1d::Int;
mass::Union{Operator,Nothing} = nothing,
θ_min::Float64 = -π / 2,
)
# setup range
dimension = operator.dimension
numbersteps = numbersteps1d^dimension
θ_max = θ_min + 2π
θ_range1d = LinRange(θ_min, θ_max, numbersteps1d)
θ_range = zeros(numbersteps, dimension)
# setup
numbermodes = size(operator.rowmodemap)[1]
eigenvalues = zeros(ComplexF64, numbersteps, numbermodes)
eigenvectors = zeros(ComplexF64, numbersteps, numbermodes, numbermodes)
rangeiterator = []
if dimension == 1
rangeiterator = θ_range1d
elseif dimension == 2
rangeiterator = Iterators.product(θ_range1d, θ_range1d)
elseif dimension == 3
rangeiterator = Iterators.product(θ_range1d, θ_range1d, θ_range1d)
end
# compute
for (step, θ) in enumerate(rangeiterator)
θ = collect(θ)
A = computesymbols(operator, θ)
if mass != nothing
A = computesymbols(mass, θ) \ A
end
currenteigen = eigen(A)
eigenvalues[step, :] = currenteigen.values
eigenvectors[step, :, :] = currenteigen.vectors
θ_range[step, :] = θ
end
# return
return (θ_range, eigenvalues, eigenvectors)
end
"""
```julia
computesymbolsoverrange(preconditioner, ω, θ; mass=nothing, θ_min=-π/2)
```
Compute the eigenvalues and eigenvectors of the symbol matrix for a preconditioned
operator over a range of θ from -π/2 to 3π/2
# Arguments:
- `preconditioner`: Preconditioner to compute symbol matries for
- `ω`: Smoothing parameter array
- `numbersteps1d`: Number of values of θ to sample in each dimension
Note: numbersteps1d^dimension symbol matrices will be computed
# Keyword Arguments:
- `mass`: Mass operator to invert for comparison to analytic solution,
default: nothing
- `θ_min`: Bottom of range of θ, shifts range to [θ_min, θ_min + 2π],
default: -π / 2
# Returns:
- Values of θ sampled
- Eigenvalues of symbol matrix at θ sampled
- Eigenvectors of symbol matrix at θ sampled
# Example:
```jldoctest
using LFAToolkit;
numbersteps1d = 5
for dimension in 1:3
# setup
mesh = []
if dimension == 1
mesh = Mesh1D(1.0);
elseif dimension == 2
mesh = Mesh2D(1.0, 1.0);
elseif dimension == 3
mesh = Mesh3D(1.0, 1.0, 1.0);
end
diffusion = GalleryOperator("diffusion", 3, 3, mesh);
# preconditioner
chebyshev = Chebyshev(diffusion);
# compute symbols
(_, eigenvalues, _) = computesymbolsoverrange(chebyshev, [1], numbersteps1d);
# verify
eigenvalues = real(eigenvalues);
if dimension == 1
@assert minimum(eigenvalues[4, :]) ≈ 0.15151515151515105
@assert maximum(eigenvalues[4, :]) ≈ 0.27272727272727226
elseif dimension == 2
@assert minimum(eigenvalues[19, :]) ≈ -0.25495098334134725
@assert maximum(eigenvalues[19, :]) ≈ -0.17128758445192374
elseif dimension == 3
@assert minimum(eigenvalues[94, :]) ≈ -0.8181818181818181
@assert maximum(eigenvalues[94, :]) ≈ -0.357575757575757
end
end
# output
```
"""
function computesymbolsoverrange(
preconditioner::AbstractPreconditioner,
ω::Array,
numbersteps1d::Int;
mass::Union{Operator,Nothing} = nothing,
θ_min::Float64 = -π / 2,
)
# setup range
dimension = preconditioner.operator.dimension
numbersteps = numbersteps1d^dimension
θ_max = θ_min + 2π
θ_range1d = LinRange(θ_min, θ_max, numbersteps1d)
θ_range = zeros(numbersteps, dimension)
# setup
numbermodes = size(preconditioner.operator.rowmodemap)[1]
eigenvalues = zeros(ComplexF64, numbersteps, numbermodes)
eigenvectors = zeros(ComplexF64, numbersteps, numbermodes, numbermodes)
rangeiterator = []
if dimension == 1
rangeiterator = θ_range1d
elseif dimension == 2
rangeiterator = Iterators.product(θ_range1d, θ_range1d)
elseif dimension == 3
rangeiterator = Iterators.product(θ_range1d, θ_range1d, θ_range1d)
end
# compute
for (step, θ) in enumerate(rangeiterator)
θ = collect(θ)
A = computesymbols(preconditioner, ω, θ)
if mass != nothing
A = computesymbols(mass, θ) \ A
end
currenteigen = eigen(A)
eigenvalues[step, :] = currenteigen.values
eigenvectors[step, :, :] = currenteigen.vectors
θ_range[step, :] = θ
end
# return
return (θ_range, eigenvalues, eigenvectors)
end
"""
```julia
computesymbolsoverrange(multigrid, ω, v, θ; mass=nothing, θ_min=-π/2)
```
Compute the eigenvalues and eigenvectors of the symbol matrix for a multigrid
preconditioned operator over a range of θ from -π/2 to 3π/2
# Arguments:
- `multigrid`: Preconditioner to compute symbol matries for
- `p`: Smoothing parameter array
- `v`: Pre and post smooths iteration count array, 0 indicates no pre or post smoothing
- `numbersteps1d`: Number of values of θ to sample in each dimension
Note: numbersteps1d^dimension symbol matrices will be computed
# Keyword Arguments:
- `mass`: Mass operator to invert for comparison to analytic solution,
default: nothing
- `θ_min`: Bottom of range of θ, shifts range to [θ_min, θ_min + 2π],
default: -π / 2
# Returns:
- Values of θ sampled
- Eigenvalues of symbol matrix at θ sampled
- Eigenvectors of symbol matrix at θ sampled
# Example:
```jldoctest
using LFAToolkit;
numbersteps1d = 5
using LinearAlgebra
for dimension in 1:3
# setup
mesh = []
if dimension == 1
mesh = Mesh1D(1.0);
elseif dimension == 2
mesh = Mesh2D(1.0, 1.0);
elseif dimension == 3
mesh = Mesh3D(1.0, 1.0, 1.0);
end
ctofbasis = TensorH1LagrangeBasis(3, 5, 1, dimension; collocatedquadrature=true);
# operators
finediffusion = GalleryOperator("diffusion", 5, 5, mesh);
coarsediffusion = GalleryOperator("diffusion", 3, 5, mesh);
# smoother
jacobi = Jacobi(finediffusion);
# preconditioner
multigrid = PMultigrid(finediffusion, coarsediffusion, jacobi, [ctofbasis]);
# compute symbols
(_, eigenvalues, _) = computesymbolsoverrange(multigrid, [1.0], [1, 1], numbersteps1d);
# verify
eigenvalues = real(eigenvalues)
if dimension == 1
@assert maximum(eigenvalues[4, :]) ≈ 0.64
elseif dimension == 2
@assert maximum(eigenvalues[19, :]) ≈ 0.9082562365654528
elseif dimension == 3
@assert maximum(eigenvalues[94, :]) ≈ 1.4359882222222669
end
end
# output
```
"""
function computesymbolsoverrange(
multigrid::Multigrid,
p::Array,
v::Array,
numbersteps1d::Int;
mass::Union{Operator,Nothing} = nothing,
θ_min::Float64 = -π / 2,
)
# setup range
dimension = multigrid.fineoperator.dimension
numbersteps = numbersteps1d^dimension
θ_max = θ_min + 2π
θ_range1d = LinRange(θ_min, θ_max, numbersteps1d)
θ_range = zeros(numbersteps, dimension)
# setup
numbermodes = size(multigrid.fineoperator.rowmodemap)[1]
eigenvalues = zeros(ComplexF64, numbersteps, numbermodes)
eigenvectors = zeros(ComplexF64, numbersteps, numbermodes, numbermodes)
rangeiterator = []
if dimension == 1
rangeiterator = θ_range1d
elseif dimension == 2
rangeiterator = Iterators.product(θ_range1d, θ_range1d)
elseif dimension == 3
rangeiterator = Iterators.product(θ_range1d, θ_range1d, θ_range1d)
end
# compute
for (step, θ) in enumerate(rangeiterator)
θ = [abs(θ_i) > 100 * eps() ? θ_i : 100 * eps() for θ_i in θ]
A = computesymbols(multigrid, p, v, θ)
if mass != nothing
A = computesymbols(mass, θ) \ A
end
currenteigen = eigen(A)
eigenvalues[step, :] = currenteigen.values
eigenvectors[step, :, :] = currenteigen.vectors
θ_range[step, :] = θ
end
# return
return (θ_range, eigenvalues, eigenvectors)
end
# ------------------------------------------------------------------------------