-
Notifications
You must be signed in to change notification settings - Fork 38
/
Copy pathksp.jl
322 lines (268 loc) · 9.9 KB
/
ksp.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
const CKSP = Ptr{Cvoid}
const CKSPType = Cstring
abstract type AbstractKSP{T, PetscLib} <: Factorization{T} end
Base.@kwdef mutable struct KSP{T, PetscLib} <: AbstractKSP{T, PetscLib}
ptr::CKSP = C_NULL
opts::Options{PetscLib}
# Stuff to keep around so that they don't get gc'ed
_A = nothing
_P = nothing
_dm = nothing
# Function pointers
ComputeRHS! = nothing
ComputeOperators! = nothing
end
struct WrappedKSP{T, PetscLib} <: AbstractKSP{T, PetscLib}
ptr::CKSP
end
scalartype(::KSP{T}) where {T} = T
Base.eltype(::KSP{T}) where {T} = T
LinearAlgebra.transpose(ksp) = LinearAlgebra.Transpose(ksp)
LinearAlgebra.adjoint(ksp) = LinearAlgebra.Adjoint(ksp)
"""
KSPSetComputeRHS!(
ksp::KSP{Number},
ComputeRHS!,
ctx = C_NULL,
)
Set the right-hand side function `ComputeRHS!` for the `ksp` using the user
`ctx`. `ComputeRHS!` should be callable with three arguments of type
`(::KSP{Number}, ::Vec, ::Ptr)`;
see [PETSc manual](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetComputeRHS.html)
"""
#function KSPSetComputeRHS! end
"""
struct Fn_KSPComputeRHS{T} end
Type used to wrap `ComputeRHS!` functions in KSP
"""
struct Fn_KSPComputeRHS{T} end
"""
KSPSetComputeOperators!(
ksp::KSP{Number},
ComputeOperators!,
ctx = C_NULL,
)
Set the linear operators function `ComputeOperators!` for the `ksp` using the
user `ctx`. `ComputeOperators!` should be callable with four arguments of type
`(::KSP{Number}, ::Mat, ::Mat, ::Ptr)`;
see [PETSc manual](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetComputeOperators.html)
"""
#function KSPSetComputeOperators! end
"""
struct Fn_KSPComputeOperators{T} end
Type used to wrap `ComputeOperators!` functions in KSP
"""
struct Fn_KSPComputeOperators{T} end
@for_libpetsc begin
function KSP{$PetscScalar}(comm::MPI.Comm; kwargs...)
@assert initialized($petsclib)
opts = Options($petsclib, kwargs...)
ksp = KSP{$PetscScalar, $PetscLib}(opts=opts)
with(ksp.opts) do
@chk ccall((:KSPCreate, $libpetsc), PetscErrorCode, (MPI.MPI_Comm, Ptr{CKSP}), comm, ksp)
end
if comm == MPI.COMM_SELF
finalizer(destroy, ksp)
end
return ksp
end
function destroy(ksp::KSP{$PetscScalar})
finalized($petsclib) ||
@chk ccall((:KSPDestroy, $libpetsc), PetscErrorCode, (Ptr{CKSP},), ksp)
return nothing
end
function setoperators!(ksp::KSP{$PetscScalar}, A::AbstractMat{$PetscScalar}, P::AbstractMat{$PetscScalar})
@chk ccall((:KSPSetOperators, $libpetsc), PetscErrorCode, (CKSP, CMat, CMat), ksp, A, P)
ksp._A = A
ksp._P = P
return nothing
end
function (::Fn_KSPComputeRHS{$PetscScalar})(
new_ksp_ptr::CKSP,
cb::CVec,
ksp_ptr::Ptr{Cvoid}
)::$PetscInt
ksp = unsafe_pointer_to_objref(ksp_ptr)
new_ksp = WrappedKSP{$PetscScalar, $PetscLib}(new_ksp_ptr)
b = Vec{$PetscScalar}(cb)
ierr = ksp.ComputeRHS!(new_ksp, b)
return $PetscInt(ierr)
end
function KSPSetComputeRHS!(
ksp::KSP{$PetscScalar},
ComputeRHS!
)
ptr_ksp = pointer_from_objref(ksp)
fptr = @cfunction(Fn_KSPComputeRHS{$PetscScalar}(),
$PetscInt,
(CKSP, CVec, Ptr{Cvoid}))
with(ksp.opts) do
@chk ccall((:KSPSetComputeRHS, $libpetsc), PetscErrorCode,
(CKSP, Ptr{Cvoid}, Ptr{Cvoid}),
ksp, fptr, ptr_ksp)
end
ksp.ComputeRHS! = ComputeRHS!
return ksp
end
function (::Fn_KSPComputeOperators{$PetscScalar})(
new_ksp_ptr::CKSP,
cA::CMat,
cP::CMat,
ksp_ptr::Ptr{Cvoid}
)::$PetscInt
ksp = unsafe_pointer_to_objref(ksp_ptr)
new_ksp = WrappedKSP{$PetscScalar, $PetscLib}(new_ksp_ptr)
A = Mat{$PetscScalar}(cA)
P = Mat{$PetscScalar}(cP)
ierr = ksp.ComputeOperators!(new_ksp, A, P)
return $PetscInt(ierr)
end
function KSPSetComputeOperators!(
ksp::KSP{$PetscScalar},
ComputeOperators!
)
ptr_ksp = pointer_from_objref(ksp)
fptr = @cfunction(Fn_KSPComputeOperators{$PetscScalar}(),
$PetscInt,
(CKSP, CMat, CMat, Ptr{Cvoid}))
with(ksp.opts) do
@chk ccall((:KSPSetComputeOperators, $libpetsc), PetscErrorCode,
(CKSP, Ptr{Cvoid}, Ptr{Cvoid}),
ksp, fptr, ptr_ksp)
end
ksp.ComputeOperators! = ComputeOperators!
return ksp
end
function KSPSetDM!(ksp::KSP{$PetscScalar}, dm::AbstractDM{$PetscLib})
with(ksp.opts) do
@chk ccall((:KSPSetDM, $libpetsc), PetscErrorCode, (CKSP, CDM), ksp, dm)
end
ksp._dm = dm
return nothing
end
function DMDA(ksp::AbstractKSP{$PetscScalar})
t_dm = Ref{CDM}()
@chk ccall(
(:KSPGetDM, $libpetsc),
PetscErrorCode,
(CKSP, Ptr{CDM}),
ksp,
t_dm,
)
dm = DMDA{$PetscLib}(t_dm[])
return dm
end
function settolerances!(ksp::KSP{$PetscScalar}; rtol=PETSC_DEFAULT, atol=PETSC_DEFAULT, divtol=PETSC_DEFAULT, max_it=PETSC_DEFAULT)
@chk ccall((:KSPSetTolerances, $libpetsc), PetscErrorCode,
(CKSP, $PetscReal, $PetscReal, $PetscReal, $PetscInt),
ksp, rtol, atol, divtol, max_it)
return nothing
end
function setfromoptions!(ksp::KSP{$PetscScalar})
with(ksp.opts) do
@chk ccall((:KSPSetFromOptions, $libpetsc), PetscErrorCode, (CKSP,), ksp)
end
end
function gettype(ksp::KSP{$PetscScalar})
t_r = Ref{CKSPType}()
@chk ccall((:KSPGetType, $libpetsc), PetscErrorCode, (CKSP, Ptr{CKSPType}), ksp, t_r)
return unsafe_string(t_r[])
end
function iters(ksp::KSP{$PetscScalar})
r_its = Ref{$PetscInt}()
@chk ccall((:KSPGetIterationNumber, $libpetsc), PetscErrorCode,
(CKSP, Ptr{$PetscInt}), ksp, r_its)
return r_its[]
end
function view(ksp::KSP{$PetscScalar}, viewer::AbstractViewer{$PetscLib}=ViewerStdout($petsclib, getcomm(ksp)))
@chk ccall((:KSPView, $libpetsc), PetscErrorCode,
(CKSP, CPetscViewer),
ksp, viewer);
return nothing
end
function resnorm(ksp::KSP{$PetscScalar})
r_rnorm = Ref{$PetscReal}()
@chk ccall((:KSPGetResidualNorm, $libpetsc), PetscErrorCode,
(CKSP, Ptr{$PetscReal}), ksp, r_rnorm)
return r_rnorm[]
end
function solve!(x::AbstractVec{$PetscScalar}, ksp::KSP{$PetscScalar, LT}, b::AbstractVec{$PetscScalar}) where LT
with(ksp.opts) do
@chk ccall((:KSPSolve, $libpetsc), PetscErrorCode,
(CKSP, CVec, CVec), ksp, b, x)
end
return x
end
function solve!(ksp::KSP{$PetscScalar})
with(ksp.opts) do
@chk ccall((:KSPSolve, $libpetsc), PetscErrorCode,
(CKSP, CVec, CVec), ksp, C_NULL, C_NULL)
end
return nothing
end
function solve!(x::AbstractVec{$PetscScalar}, tksp::Transpose{T,K}, b::AbstractVec{$PetscScalar}) where {T,K <: KSP{$PetscScalar}}
ksp = parent(tksp)
with(ksp.opts) do
@chk ccall((:KSPSolveTranspose, $libpetsc), PetscErrorCode,
(CKSP, CVec, CVec), ksp, b, x)
end
return x
end
function solve!(x::AbstractVec{$PetscScalar}, tksp::LinearAlgebra.AdjointFactorization{T,K}, b::AbstractVec{$PetscScalar}) where {T,K <: KSP{$PetscScalar}}
ksp = parent(tksp)
with(ksp.opts) do
@chk ccall((:KSPSolveTranspose, $libpetsc), PetscErrorCode,
(CKSP, CVec, CVec), ksp, b, x)
end
return x
end
end
# no generic Adjoint solve defined, but for Real we can use Adjoint
solve!(x::AbstractVec{T}, aksp::Adjoint{T,K}, b::AbstractVec{T}) where {K <: KSP{T}} where {T<:Real} =
solve!(x, transpose(parent(aksp)), b)
const KSPAT{T, LT} = Union{KSP{T, LT}, Transpose{T, KSP{T, LT}}, Adjoint{T, KSP{T, LT}}, LinearAlgebra.AdjointFactorization{T, KSP{T, LT}}}
LinearAlgebra.ldiv!(x::AbstractVec{T}, ksp::KSPAT{T, LT}, b::AbstractVec{T}) where {T, LT} = solve!(x, ksp, b)
function LinearAlgebra.ldiv!(x::AbstractVector{T}, ksp::KSPAT{T, LT}, b::AbstractVector{T}) where {T, LT}
parent(solve!(AbstractVec(x), ksp, AbstractVec(b)))
end
Base.:\(ksp::KSPAT{T, LT}, b::AbstractVector{T}) where {T, LT} = ldiv!(similar(b), ksp, b)
#Base.:\(kspt::LinearAlgebra.AdjointFactorization{T,KSPT{T, LT}}, b::AbstractVector{T}) where {T, LT} = ldiv!(similar(b), kspt._A', b)
"""
KSP(A, P; options...)
Construct a PETSc Krylov subspace solver.
Any PETSc options prefixed with `ksp_` and `pc_` can be passed as keywords.
"""
function KSP(A::AbstractMat{T}, P::AbstractMat{T}=A; kwargs...) where {T}
ksp = KSP{T}(getcomm(A); kwargs...)
setoperators!(ksp, A, P)
setfromoptions!(ksp)
return ksp
end
"""
KSP(da::AbstractDM; options...)
Construct a PETSc Krylov subspace solver from the distributed mesh
Any PETSc options prefixed with `ksp_` and `pc_` can be passed as keywords.
see [PETSc manual](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetDM.html)
"""
function KSP(dm::AbstractDM{PetscLib}; kwargs...) where {PetscLib}
T = scalartype(PetscLib)
ksp = KSP{T}(getcomm(dm); kwargs...)
KSPSetDM!(ksp, dm)
setfromoptions!(ksp)
return ksp
end
Base.show(io::IO, ksp::KSP) = _show(io, ksp)
"""
iters(ksp::KSP)
Gets the current iteration number; if the `solve!` is complete, returns the number of iterations used.
# External Links
$(_doc_external("KSP/KSPGetIterationNumber"))
"""
iters
"""
resnorm(ksp::KSP)
Gets the last (approximate preconditioned) residual norm that has been computed.
# External Links
$(_doc_external("KSP/KSPGetResidualNorm"))
"""
resnorm