-
Notifications
You must be signed in to change notification settings - Fork 195
/
conjugate_gradient_poisson_solver.jl
177 lines (141 loc) · 7.59 KB
/
conjugate_gradient_poisson_solver.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
using Oceananigans.Operators: divᶜᶜᶜ, ∇²ᶜᶜᶜ
using Statistics: mean
using KernelAbstractions: @kernel, @index
import Oceananigans.Architectures: architecture
struct ConjugateGradientPoissonSolver{G, R, S}
grid :: G
right_hand_side :: R
conjugate_gradient_solver :: S
end
architecture(solver::ConjugateGradientPoissonSolver) = architecture(cgps.grid)
iteration(cgps::ConjugateGradientPoissonSolver) = iteration(cgps.conjugate_gradient_solver)
Base.summary(ips::ConjugateGradientPoissonSolver) =
summary("ConjugateGradientPoissonSolver on ", summary(ips.grid))
function Base.show(io::IO, ips::ConjugateGradientPoissonSolver)
A = architecture(ips.grid)
print(io, "ConjugateGradientPoissonSolver:", '\n',
"├── grid: ", summary(ips.grid), '\n',
"└── conjugate_gradient_solver: ", summary(ips.conjugate_gradient_solver), '\n',
" ├── maxiter: ", prettysummary(ips.conjugate_gradient_solver.maxiter), '\n',
" ├── reltol: ", prettysummary(ips.conjugate_gradient_solver.reltol), '\n',
" ├── abstol: ", prettysummary(ips.conjugate_gradient_solver.abstol), '\n',
" ├── preconditioner: ", prettysummary(ips.conjugate_gradient_solver.preconditioner), '\n',
" └── iteration: ", prettysummary(ips.conjugate_gradient_solver.iteration))
end
@kernel function laplacian!(∇²ϕ, grid, ϕ)
i, j, k = @index(Global, NTuple)
@inbounds ∇²ϕ[i, j, k] = ∇²ᶜᶜᶜ(i, j, k, grid, ϕ)
end
function compute_laplacian!(∇²ϕ, ϕ)
grid = ϕ.grid
arch = architecture(grid)
fill_halo_regions!(ϕ)
launch!(arch, grid, :xyz, laplacian!, ∇²ϕ, grid, ϕ)
return nothing
end
struct DefaultPreconditioner end
function ConjugateGradientPoissonSolver(grid;
preconditioner = DefaultPreconditioner(),
reltol = sqrt(eps(grid)),
abstol = sqrt(eps(grid)),
kw...)
if preconditioner isa DefaultPreconditioner # try to make a useful default
if grid isa ImmersedBoundaryGrid && grid.underlying_grid isa GridWithFFTSolver
preconditioner = fft_poisson_solver(grid.underlying_grid)
else
preconditioner = DiagonallyDominantPreconditioner()
end
end
rhs = CenterField(grid)
conjugate_gradient_solver = ConjugateGradientSolver(compute_laplacian!;
reltol,
abstol,
preconditioner,
template_field = rhs,
kw...)
return ConjugateGradientPoissonSolver(grid, rhs, conjugate_gradient_solver)
end
#####
##### A preconditioner based on the FFT solver
#####
@kernel function fft_preconditioner_rhs!(preconditioner_rhs, rhs)
i, j, k = @index(Global, NTuple)
@inbounds preconditioner_rhs[i, j, k] = rhs[i, j, k]
end
@kernel function fourier_tridiagonal_preconditioner_rhs!(preconditioner_rhs, ::XDirection, grid, rhs)
i, j, k = @index(Global, NTuple)
@inbounds preconditioner_rhs[i, j, k] = Δxᶜᶜᶜ(i, j, k, grid) * rhs[i, j, k]
end
@kernel function fourier_tridiagonal_preconditioner_rhs!(preconditioner_rhs, ::YDirection, grid, rhs)
i, j, k = @index(Global, NTuple)
@inbounds preconditioner_rhs[i, j, k] = Δyᶜᶜᶜ(i, j, k, grid) * rhs[i, j, k]
end
@kernel function fourier_tridiagonal_preconditioner_rhs!(preconditioner_rhs, ::ZDirection, grid, rhs)
i, j, k = @index(Global, NTuple)
@inbounds preconditioner_rhs[i, j, k] = Δzᶜᶜᶜ(i, j, k, grid) * rhs[i, j, k]
end
function compute_preconditioner_rhs!(solver::FFTBasedPoissonSolver, rhs)
grid = solver.grid
arch = architecture(grid)
launch!(arch, grid, :xyz, fft_preconditioner_rhs!, solver.storage, rhs)
return nothing
end
function compute_preconditioner_rhs!(solver::FourierTridiagonalPoissonSolver, rhs)
grid = solver.grid
arch = architecture(grid)
tridiagonal_dir = solver.batched_tridiagonal_solver.tridiagonal_direction
launch!(arch, grid, :xyz, fourier_tridiagonal_preconditioner_rhs!,
solver.storage, tridiagonal_dir, rhs)
return nothing
end
const FFTBasedPreconditioner = Union{FFTBasedPoissonSolver, FourierTridiagonalPoissonSolver}
function precondition!(p, preconditioner::FFTBasedPreconditioner, r, args...)
compute_preconditioner_rhs!(preconditioner, r)
p = solve!(p, preconditioner)
mean_p = mean(p)
grid = p.grid
arch = architecture(grid)
launch!(arch, grid, :xyz, subtract_and_mask!, p, grid, mean_p)
return p
end
@kernel function subtract_and_mask!(a, grid, b)
i, j, k = @index(Global, NTuple)
active = !inactive_cell(i, j, k, grid)
a[i, j, k] = (a[i, j, k] - b) * active
end
#####
##### The "DiagonallyDominantPreconditioner" (Marshall et al 1997)
#####
struct DiagonallyDominantPreconditioner end
Base.summary(::DiagonallyDominantPreconditioner) = "DiagonallyDominantPreconditioner"
@inline function precondition!(p, ::DiagonallyDominantPreconditioner, r, args...)
grid = r.grid
arch = architecture(p)
fill_halo_regions!(r)
launch!(arch, grid, :xyz, _diagonally_dominant_precondition!, p, grid, r)
mean_p = mean(p)
launch!(arch, grid, :xyz, subtract_and_mask!, p, grid, mean_p)
return p
end
# Kernels that calculate coefficients for the preconditioner
@inline Ax⁻(i, j, k, grid) = Axᶠᶜᶜ(i, j, k, grid) / Δxᶠᶜᶜ(i, j, k, grid) / Vᶜᶜᶜ(i, j, k, grid)
@inline Ax⁺(i, j, k, grid) = Axᶠᶜᶜ(i+1, j, k, grid) / Δxᶠᶜᶜ(i+1, j, k, grid) / Vᶜᶜᶜ(i, j, k, grid)
@inline Ay⁻(i, j, k, grid) = Ayᶜᶠᶜ(i, j, k, grid) / Δyᶜᶠᶜ(i, j, k, grid) / Vᶜᶜᶜ(i, j, k, grid)
@inline Ay⁺(i, j, k, grid) = Ayᶜᶠᶜ(i, j+1, k, grid) / Δyᶜᶠᶜ(i, j+1, k, grid) / Vᶜᶜᶜ(i, j, k, grid)
@inline Az⁻(i, j, k, grid) = Azᶜᶜᶠ(i, j, k, grid) / Δzᶜᶜᶠ(i, j, k, grid) / Vᶜᶜᶜ(i, j, k, grid)
@inline Az⁺(i, j, k, grid) = Azᶜᶜᶠ(i, j, k+1, grid) / Δzᶜᶜᶠ(i, j, k+1, grid) / Vᶜᶜᶜ(i, j, k, grid)
@inline Ac(i, j, k, grid) = - Ax⁻(i, j, k, grid) - Ax⁺(i, j, k, grid) -
Ay⁻(i, j, k, grid) - Ay⁺(i, j, k, grid) -
Az⁻(i, j, k, grid) - Az⁺(i, j, k, grid)
@inline heuristic_residual(i, j, k, grid, r) =
@inbounds 1 / Ac(i, j, k, grid) * (r[i, j, k] - 2 * Ax⁻(i, j, k, grid) / (Ac(i, j, k, grid) + Ac(i-1, j, k, grid)) * r[i-1, j, k] -
2 * Ax⁺(i, j, k, grid) / (Ac(i, j, k, grid) + Ac(i+1, j, k, grid)) * r[i+1, j, k] -
2 * Ay⁻(i, j, k, grid) / (Ac(i, j, k, grid) + Ac(i, j-1, k, grid)) * r[i, j-1, k] -
2 * Ay⁺(i, j, k, grid) / (Ac(i, j, k, grid) + Ac(i, j+1, k, grid)) * r[i, j+1, k] -
2 * Az⁻(i, j, k, grid) / (Ac(i, j, k, grid) + Ac(i, j, k-1, grid)) * r[i, j, k-1] -
2 * Az⁺(i, j, k, grid) / (Ac(i, j, k, grid) + Ac(i, j, k+1, grid)) * r[i, j, k+1])
@kernel function _diagonally_dominant_precondition!(p, grid, r)
i, j, k = @index(Global, NTuple)
active = !inactive_cell(i, j, k, grid)
@inbounds p[i, j, k] = heuristic_residual(i, j, k, grid, r) * active
end