-
Notifications
You must be signed in to change notification settings - Fork 209
/
preconditioners.jl
328 lines (305 loc) · 15.1 KB
/
preconditioners.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
# routines that implement different preconditioners
export ic02!, ic02, ilu02!, ilu02, gtsv2!, gtsv2
"""
ic02!(A::CuSparseMatrix, index::SparseChar='O')
Incomplete Cholesky factorization with no pivoting.
Preserves the sparse layout of matrix `A`.
"""
function ic02! end
"""
ilu02!(A::CuSparseMatrix, index::SparseChar='O')
Incomplete LU factorization with no pivoting.
Preserves the sparse layout of matrix `A`.
"""
function ilu02! end
"""
gtsv2!(dl::CuVector, d::CuVector, du::CuVector, B::CuVecOrMat, index::SparseChar='O'; pivoting::Bool=true)
Solve the linear system `A * X = B` where `A` is a tridiagonal matrix defined
by three vectors corresponding to its lower (`dl`), main (`d`), and upper (`du`) diagonals.
With `pivoting`, the solution is more accurate but also more expensive.
Note that the solution `X` overwrites the right-hand side `B`.
"""
function gtsv2! end
# csric02
for (bname,aname,sname,elty) in ((:cusparseScsric02_bufferSize, :cusparseScsric02_analysis, :cusparseScsric02, :Float32),
(:cusparseDcsric02_bufferSize, :cusparseDcsric02_analysis, :cusparseDcsric02, :Float64),
(:cusparseCcsric02_bufferSize, :cusparseCcsric02_analysis, :cusparseCcsric02, :ComplexF32),
(:cusparseZcsric02_bufferSize, :cusparseZcsric02_analysis, :cusparseZcsric02, :ComplexF64))
@eval begin
function ic02!(A::CuSparseMatrixCSR{$elty}, index::SparseChar='O')
desc = CuMatrixDescriptor('G', 'L', 'N', index)
m,n = size(A)
if m != n
throw(DimensionMismatch("A must be square, but has dimensions ($m,$n)!"))
end
info = IC0Info()
function bufferSize()
out = Ref{Cint}(1)
$bname(handle(), m, nnz(A), desc, nonzeros(A), A.rowPtr, A.colVal, info,
out)
return out[]
end
with_workspace(bufferSize) do buffer
$aname(handle(), m, nnz(A), desc,
nonzeros(A), A.rowPtr, A.colVal, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
posit = Ref{Cint}(1)
cusparseXcsric02_zeroPivot(handle(), info, posit)
if posit[] >= 0
error("Structural/numerical zero in A at ($(posit[]),$(posit[])))")
end
$sname(handle(), m, nnz(A),
desc, nonzeros(A), A.rowPtr, A.colVal, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
end
A
end
end
end
# cscic02
for (bname,aname,sname,elty) in ((:cusparseScsric02_bufferSize, :cusparseScsric02_analysis, :cusparseScsric02, :Float32),
(:cusparseDcsric02_bufferSize, :cusparseDcsric02_analysis, :cusparseDcsric02, :Float64),
(:cusparseCcsric02_bufferSize, :cusparseCcsric02_analysis, :cusparseCcsric02, :ComplexF32),
(:cusparseZcsric02_bufferSize, :cusparseZcsric02_analysis, :cusparseZcsric02, :ComplexF64))
@eval begin
function ic02!(A::CuSparseMatrixCSC{$elty}, index::SparseChar='O')
desc = CuMatrixDescriptor('G', 'L', 'N', index)
m,n = size(A)
if m != n
throw(DimensionMismatch("A must be square, but has dimensions ($m,$n)!"))
end
info = IC0Info()
function bufferSize()
out = Ref{Cint}(1)
$bname(handle(), m, nnz(A), desc, nonzeros(A), A.colPtr, rowvals(A),
info, out)
return out[]
end
with_workspace(bufferSize) do buffer
$aname(handle(), m, nnz(A), desc,
nonzeros(A), A.colPtr, rowvals(A), info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
posit = Ref{Cint}(1)
cusparseXcsric02_zeroPivot(handle(), info, posit)
if posit[] >= 0
error("Structural/numerical zero in A at ($(posit[]),$(posit[])))")
end
$sname(handle(), m, nnz(A),
desc, nonzeros(A), A.colPtr, rowvals(A), info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
end
A
end
end
end
# csrilu02
for (bname,aname,sname,elty) in ((:cusparseScsrilu02_bufferSize, :cusparseScsrilu02_analysis, :cusparseScsrilu02, :Float32),
(:cusparseDcsrilu02_bufferSize, :cusparseDcsrilu02_analysis, :cusparseDcsrilu02, :Float64),
(:cusparseCcsrilu02_bufferSize, :cusparseCcsrilu02_analysis, :cusparseCcsrilu02, :ComplexF32),
(:cusparseZcsrilu02_bufferSize, :cusparseZcsrilu02_analysis, :cusparseZcsrilu02, :ComplexF64))
@eval begin
function ilu02!(A::CuSparseMatrixCSR{$elty}, index::SparseChar='O')
desc = CuMatrixDescriptor('G', 'L', 'N', index)
m,n = size(A)
if m != n
throw(DimensionMismatch("A must be square, but has dimensions ($m,$n)!"))
end
info = ILU0Info()
function bufferSize()
out = Ref{Cint}(1)
$bname(handle(), m, nnz(A), desc,
nonzeros(A), A.rowPtr, A.colVal, info,
out)
return out[]
end
with_workspace(bufferSize) do buffer
$aname(handle(), m, nnz(A), desc,
nonzeros(A), A.rowPtr, A.colVal, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
posit = Ref{Cint}(1)
cusparseXcsrilu02_zeroPivot(handle(), info, posit)
if posit[] >= 0
error("Structural zero in A at ($(posit[]),$(posit[])))")
end
$sname(handle(), m, nnz(A),
desc, nonzeros(A), A.rowPtr, A.colVal, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
end
A
end
end
end
# cscilu02
for (bname,aname,sname,elty) in ((:cusparseScsrilu02_bufferSize, :cusparseScsrilu02_analysis, :cusparseScsrilu02, :Float32),
(:cusparseDcsrilu02_bufferSize, :cusparseDcsrilu02_analysis, :cusparseDcsrilu02, :Float64),
(:cusparseCcsrilu02_bufferSize, :cusparseCcsrilu02_analysis, :cusparseCcsrilu02, :ComplexF32),
(:cusparseZcsrilu02_bufferSize, :cusparseZcsrilu02_analysis, :cusparseZcsrilu02, :ComplexF64))
@eval begin
function ilu02!(A::CuSparseMatrixCSC{$elty}, index::SparseChar='O')
desc = CuMatrixDescriptor('G', 'L', 'N', index)
m,n = size(A)
if m != n
throw(DimensionMismatch("A must be square, but has dimensions ($m,$n)!"))
end
info = ILU0Info()
function bufferSize()
out = Ref{Cint}(1)
$bname(handle(), m, nnz(A), desc,
nonzeros(A), A.colPtr, rowvals(A), info,
out)
return out[]
end
with_workspace(bufferSize) do buffer
$aname(handle(), m, nnz(A), desc,
nonzeros(A), A.colPtr, rowvals(A), info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
posit = Ref{Cint}(1)
cusparseXcsrilu02_zeroPivot(handle(), info, posit)
if posit[] >= 0
error("Structural zero in A at ($(posit[]),$(posit[])))")
end
$sname(handle(), m, nnz(A),
desc, nonzeros(A), A.colPtr, rowvals(A), info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
end
A
end
end
end
# bsric02
for (bname,aname,sname,elty) in ((:cusparseSbsric02_bufferSize, :cusparseSbsric02_analysis, :cusparseSbsric02, :Float32),
(:cusparseDbsric02_bufferSize, :cusparseDbsric02_analysis, :cusparseDbsric02, :Float64),
(:cusparseCbsric02_bufferSize, :cusparseCbsric02_analysis, :cusparseCbsric02, :ComplexF32),
(:cusparseZbsric02_bufferSize, :cusparseZbsric02_analysis, :cusparseZbsric02, :ComplexF64))
@eval begin
function ic02!(A::CuSparseMatrixBSR{$elty}, index::SparseChar='O')
desc = CuMatrixDescriptor('G', 'U', 'N', index)
m,n = size(A)
if m != n
throw(DimensionMismatch("A must be square, but has dimensions ($m,$n)!"))
end
mb = div(m,A.blockDim)
info = IC0InfoBSR()
function bufferSize()
out = Ref{Cint}(1)
$bname(handle(), A.dir, mb, nnz(A), desc, nonzeros(A),
A.rowPtr, A.colVal, A.blockDim, info, out)
return out[]
end
with_workspace(bufferSize) do buffer
$aname(handle(), A.dir, mb, nnz(A), desc,
nonzeros(A), A.rowPtr, A.colVal, A.blockDim, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
posit = Ref{Cint}(1)
cusparseXbsric02_zeroPivot(handle(), info, posit)
if posit[] >= 0
error("Structural/numerical zero in A at ($(posit[]),$(posit[])))")
end
$sname(handle(), A.dir, mb, nnz(A), desc,
nonzeros(A), A.rowPtr, A.colVal, A.blockDim, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
end
A
end
end
end
# bsrilu02
for (bname,aname,sname,elty) in ((:cusparseSbsrilu02_bufferSize, :cusparseSbsrilu02_analysis, :cusparseSbsrilu02, :Float32),
(:cusparseDbsrilu02_bufferSize, :cusparseDbsrilu02_analysis, :cusparseDbsrilu02, :Float64),
(:cusparseCbsrilu02_bufferSize, :cusparseCbsrilu02_analysis, :cusparseCbsrilu02, :ComplexF32),
(:cusparseZbsrilu02_bufferSize, :cusparseZbsrilu02_analysis, :cusparseZbsrilu02, :ComplexF64))
@eval begin
function ilu02!(A::CuSparseMatrixBSR{$elty}, index::SparseChar='O')
desc = CuMatrixDescriptor('G', 'U', 'N', index)
m,n = size(A)
if m != n
throw(DimensionMismatch("A must be square, but has dimensions ($m,$n)!"))
end
mb = div(m,A.blockDim)
info = ILU0InfoBSR()
function bufferSize()
out = Ref{Cint}(1)
$bname(handle(), A.dir, mb, nnz(A), desc, nonzeros(A),
A.rowPtr, A.colVal, A.blockDim, info, out)
return out[]
end
with_workspace(bufferSize) do buffer
$aname(handle(), A.dir, mb, nnz(A), desc,
nonzeros(A), A.rowPtr, A.colVal, A.blockDim, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
posit = Ref{Cint}(1)
cusparseXbsrilu02_zeroPivot(handle(), info, posit)
if posit[] >= 0
error("Structural/numerical zero in A at ($(posit[]),$(posit[])))")
end
$sname(handle(), A.dir, mb, nnz(A), desc,
nonzeros(A), A.rowPtr, A.colVal, A.blockDim, info,
CUSPARSE_SOLVE_POLICY_USE_LEVEL, buffer)
end
A
end
end
end
for elty in (:Float32, :Float64, :ComplexF32, :ComplexF64)
@eval begin
function ilu02(A::CuSparseMatrix{$elty}, index::SparseChar='O')
isa(A, CuSparseMatrixCOO) && throw(ErrorException("ILU(0) decomposition of CuSparseMatrixCOO is not supported by the current CUDA version."))
ilu02!(copy(A),index)
end
function ic02(A::CuSparseMatrix{$elty}, index::SparseChar='O')
isa(A, CuSparseMatrixCOO) && throw(ErrorException("IC(0) decomposition of CuSparseMatrixCOO is not supported by the current CUDA version."))
ic02!(copy(A),index)
end
function ilu02(A::HermOrSym{$elty,CuSparseMatrix{$elty}}, index::SparseChar='O')
isa(A, CuSparseMatrixCOO) && throw(ErrorException("ILU(0) decomposition of CuSparseMatrixCOO is not supported by the current CUDA version."))
ilu02!(copy(A.data),index)
end
function ic02(A::HermOrSym{$elty,CuSparseMatrix{$elty}}, index::SparseChar='O')
isa(A, CuSparseMatrixCOO) && throw(ErrorException("IC(0) decomposition of CuSparseMatrixCOO is not supported by the current CUDA version."))
ic02!(copy(A.data),index)
end
end
end
# gtsv2
for (bname_pivot,fname_pivot,bname_nopivot,fname_nopivot,elty) in ((:cusparseSgtsv2_bufferSizeExt, :cusparseSgtsv2, :cusparseSgtsv2_nopivot_bufferSizeExt, :cusparseSgtsv2_nopivot, :Float32),
(:cusparseDgtsv2_bufferSizeExt, :cusparseDgtsv2, :cusparseDgtsv2_nopivot_bufferSizeExt, :cusparseDgtsv2_nopivot, :Float64),
(:cusparseCgtsv2_bufferSizeExt, :cusparseCgtsv2, :cusparseCgtsv2_nopivot_bufferSizeExt, :cusparseCgtsv2_nopivot, :ComplexF32),
(:cusparseZgtsv2_bufferSizeExt, :cusparseZgtsv2, :cusparseZgtsv2_nopivot_bufferSizeExt, :cusparseZgtsv2_nopivot, :ComplexF64))
@eval begin
function gtsv2!(dl::CuVector{$elty}, d::CuVector{$elty}, du::CuVector{$elty}, B::CuVecOrMat{$elty}, index::SparseChar='O'; pivoting::Bool=true)
ml = length(dl)
m = length(d)
mu = length(du)
mB = size(B,1)
(m ≤ 2) && throw(DimensionMismatch("The size of the linear system must be at least 3."))
!(ml == m == mu) && throw(DimensionMismatch("(dl, d, du) must have the same length, the size of the vectors is ($ml,$m,$mu)!"))
(m != mB) && throw(DimensionMismatch("The tridiagonal matrix and the right-hand side B have inconsistent dimensions ($m != $mB)!"))
n = size(B,2)
ldb = max(1,stride(B,2))
function bufferSize()
out = Ref{Csize_t}(1)
if pivoting
$bname_pivot(handle(), m, n, dl, d, du, B, ldb, out)
else
$bname_nopivot(handle(), m, n, dl, d, du, B, ldb, out)
end
return out[]
end
with_workspace(bufferSize) do buffer
if pivoting
$fname_pivot(handle(), m, n, dl, d, du, B, ldb, buffer)
else
$fname_nopivot(handle(), m, n, dl, d, du, B, ldb, buffer)
end
end
B
end
end
end
for elty in (:Float32, :Float64, :ComplexF32, :ComplexF64)
@eval begin
function gtsv2(dl::CuVector{$elty}, d::CuVector{$elty}, du::CuVector{$elty}, B::CuVecOrMat{$elty}, index::SparseChar='O'; pivoting::Bool=true)
gtsv2!(dl, d, du, copy(B), index; pivoting)
end
end
end