Skip to content

Commit 2c915ab

Browse files
musmandreasnoack
authored andcommitted
Update deprecated ccall syntax in arpack (#23307)
1 parent 1dedf50 commit 2c915ab

File tree

1 file changed

+92
-90
lines changed

1 file changed

+92
-90
lines changed

base/linalg/arpack.jl

Lines changed: 92 additions & 90 deletions
Original file line numberDiff line numberDiff line change
@@ -12,8 +12,7 @@ function aupd_wrapper(T, matvecA!::Function, matvecB::Function, solveSI::Functio
1212
tol::Real, maxiter::Integer, mode::Integer, v0::Vector)
1313
lworkl = cmplx ? ncv * (3*ncv + 5) : (sym ? ncv * (ncv + 8) : ncv * (3*ncv + 6) )
1414
TR = cmplx ? T.types[1] : T
15-
TOL = Vector{TR}(1)
16-
TOL[1] = tol
15+
TOL = Ref{TR}(tol)
1716

1817
v = Matrix{T}(n, ncv)
1918
workd = Vector{T}(3*n)
@@ -22,14 +21,14 @@ function aupd_wrapper(T, matvecA!::Function, matvecB::Function, solveSI::Functio
2221

2322
if isempty(v0)
2423
resid = Vector{T}(n)
25-
info = zeros(BlasInt, 1)
24+
info = Ref{BlasInt}(0)
2625
else
2726
resid = deepcopy(v0)
28-
info = ones(BlasInt, 1)
27+
info = Ref{BlasInt}(1)
2928
end
3029
iparam = zeros(BlasInt, 11)
3130
ipntr = zeros(BlasInt, (sym && !cmplx) ? 11 : 14)
32-
ido = zeros(BlasInt, 1)
31+
ido = Ref{BlasInt}(0)
3332

3433
iparam[1] = BlasInt(1) # ishifts
3534
iparam[3] = BlasInt(maxiter) # maxiter
@@ -48,50 +47,50 @@ function aupd_wrapper(T, matvecA!::Function, matvecB::Function, solveSI::Functio
4847
naupd(ido, bmat, n, which, nev, TOL, resid, ncv, v, n,
4948
iparam, ipntr, workd, workl, lworkl, info)
5049
end
51-
if info[1] != 0
52-
throw(ARPACKException(info[1]))
50+
if info[] != 0
51+
throw(ARPACKException(info[]))
5352
end
5453

5554
x = view(workd, ipntr[1]+zernm1)
5655
y = view(workd, ipntr[2]+zernm1)
5756
if mode == 1 # corresponds to dsdrv1, dndrv1 or zndrv1
58-
if ido[1] == 1
57+
if ido[] == 1
5958
matvecA!(y, x)
60-
elseif ido[1] == 99
59+
elseif ido[] == 99
6160
break
6261
else
6362
throw(ARPACKException("unexpected behavior"))
6463
end
6564
elseif mode == 3 && bmat == "I" # corresponds to dsdrv2, dndrv2 or zndrv2
66-
if ido[1] == -1 || ido[1] == 1
65+
if ido[] == -1 || ido[] == 1
6766
y[:] = solveSI(x)
68-
elseif ido[1] == 99
67+
elseif ido[] == 99
6968
break
7069
else
7170
throw(ARPACKException("unexpected behavior"))
7271
end
7372
elseif mode == 2 # corresponds to dsdrv3, dndrv3 or zndrv3
74-
if ido[1] == -1 || ido[1] == 1
73+
if ido[] == -1 || ido[] == 1
7574
matvecA!(y, x)
7675
if sym
7776
x[:] = y # overwrite as per Remark 5 in dsaupd.f
7877
end
7978
y[:] = solveSI(y)
80-
elseif ido[1] == 2
79+
elseif ido[] == 2
8180
y[:] = matvecB(x)
82-
elseif ido[1] == 99
81+
elseif ido[] == 99
8382
break
8483
else
8584
throw(ARPACKException("unexpected behavior"))
8685
end
8786
elseif mode == 3 && bmat == "G" # corresponds to dsdrv4, dndrv4 or zndrv4
88-
if ido[1] == -1
87+
if ido[] == -1
8988
y[:] = solveSI(matvecB(x))
90-
elseif ido[1] == 1
89+
elseif ido[] == 1
9190
y[:] = solveSI(view(workd,ipntr[3]+zernm1))
92-
elseif ido[1] == 2
91+
elseif ido[] == 2
9392
y[:] = matvecB(x)
94-
elseif ido[1] == 99
93+
elseif ido[] == 99
9594
break
9695
else
9796
throw(ARPACKException("unexpected behavior"))
@@ -106,63 +105,63 @@ end
106105

107106
function eupd_wrapper(T, n::Integer, sym::Bool, cmplx::Bool, bmat::String,
108107
nev::Integer, which::String, ritzvec::Bool,
109-
TOL::Array, resid, ncv::Integer, v, ldv, sigma, iparam, ipntr,
108+
TOL::Ref, resid, ncv::Integer, v, ldv, sigma, iparam, ipntr,
110109
workd, workl, lworkl, rwork)
111110
howmny = "A"
112111
select = Vector{BlasInt}(ncv)
113-
info = zeros(BlasInt, 1)
112+
info = Ref{BlasInt}(0)
114113

115-
dmap = x->abs.(x)
114+
dmap = x -> abs.(x)
116115
if iparam[7] == 3 # shift-and-invert
117-
dmap = x->abs.(1 ./ (x .- sigma))
116+
dmap = x -> abs.(1 ./ (x .- sigma))
118117
elseif which == "LR" || which == "LA" || which == "BE"
119118
dmap = real
120119
elseif which == "SR" || which == "SA"
121-
dmap = x->-real(x)
120+
dmap = x -> -real(x)
122121
elseif which == "LI"
123122
dmap = imag
124123
elseif which == "SI"
125-
dmap = x->-imag(x)
124+
dmap = x -> -imag(x)
126125
end
127126

128127
if cmplx
129128
d = Vector{T}(nev+1)
130-
sigmar = ones(T, 1)*sigma
129+
sigmar = Ref{T}(sigma)
131130
workev = Vector{T}(2ncv)
132131
neupd(ritzvec, howmny, select, d, v, ldv, sigmar, workev,
133132
bmat, n, which, nev, TOL, resid, ncv, v, ldv,
134133
iparam, ipntr, workd, workl, lworkl, rwork, info)
135-
if info[1] != 0
136-
throw(ARPACKException(info[1]))
134+
if info[] != 0
135+
throw(ARPACKException(info[]))
137136
end
138137

139138
p = sortperm(dmap(d[1:nev]), rev=true)
140139
return ritzvec ? (d[p], v[1:n, p],iparam[5],iparam[3],iparam[9],resid) : (d[p],iparam[5],iparam[3],iparam[9],resid)
141140
elseif sym
142141
d = Vector{T}(nev)
143-
sigmar = ones(T, 1)*sigma
142+
sigmar = Ref{T}(sigma)
144143
seupd(ritzvec, howmny, select, d, v, ldv, sigmar,
145144
bmat, n, which, nev, TOL, resid, ncv, v, ldv,
146145
iparam, ipntr, workd, workl, lworkl, info)
147-
if info[1] != 0
148-
throw(ARPACKException(info[1]))
146+
if info[] != 0
147+
throw(ARPACKException(info[]))
149148
end
150149

151150
p = sortperm(dmap(d), rev=true)
152151
return ritzvec ? (d[p], v[1:n, p],iparam[5],iparam[3],iparam[9],resid) : (d[p],iparam[5],iparam[3],iparam[9],resid)
153152
else
154-
dr = Vector{T}(nev+1)
155-
di = Vector{T}(nev+1)
153+
dr = Vector{T}(nev+1)
154+
di = Vector{T}(nev+1)
156155
fill!(dr,NaN)
157156
fill!(di,NaN)
158-
sigmar = ones(T, 1)*real(sigma)
159-
sigmai = ones(T, 1)*imag(sigma)
157+
sigmar = Ref{T}(real(sigma))
158+
sigmai = Ref{T}(imag(sigma))
160159
workev = Vector{T}(3*ncv)
161160
neupd(ritzvec, howmny, select, dr, di, v, ldv, sigmar, sigmai,
162161
workev, bmat, n, which, nev, TOL, resid, ncv, v, ldv,
163162
iparam, ipntr, workd, workl, lworkl, info)
164-
if info[1] != 0
165-
throw(ARPACKException(info[1]))
163+
if info[] != 0
164+
throw(ARPACKException(info[]))
166165
end
167166
evec = complex.(Matrix{T}(n, nev+1), Matrix{T}(n, nev+1))
168167

@@ -203,52 +202,54 @@ for (T, saupd_name, seupd_name, naupd_name, neupd_name) in
203202
((:Float64, :dsaupd_, :dseupd_, :dnaupd_, :dneupd_),
204203
(:Float32, :ssaupd_, :sseupd_, :snaupd_, :sneupd_))
205204
@eval begin
206-
function naupd(ido, bmat, n, evtype, nev, TOL::Array{$T}, resid::Array{$T}, ncv, v::Array{$T}, ldv,
207-
iparam, ipntr, workd::Array{$T}, workl::Array{$T}, lworkl, info)
205+
function naupd(ido, bmat, n, evtype, nev, TOL::Ref{$T}, resid::Vector{$T}, ncv, v::Matrix{$T}, ldv,
206+
iparam, ipntr, workd::Vector{$T}, workl::Vector{$T}, lworkl, info)
208207
ccall(($(string(naupd_name)), :libarpack), Void,
209-
(Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt},
210-
Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{$T}, Ptr{BlasInt},
211-
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{BlasInt}, Clong, Clong),
212-
ido, bmat, &n, evtype, &nev, TOL, resid, &ncv, v, &ldv,
213-
iparam, ipntr, workd, workl, &lworkl, info, sizeof(bmat), sizeof(evtype))
208+
(Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt},
209+
Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ptr{$T}, Ref{BlasInt},
210+
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ref{BlasInt}, Clong, Clong),
211+
ido, bmat, n, evtype, nev,
212+
TOL, resid, ncv, v, ldv,
213+
iparam, ipntr, workd, workl, lworkl, info, sizeof(bmat), sizeof(evtype))
214214
end
215215

216216
function neupd(rvec, howmny, select, dr, di, z, ldz, sigmar, sigmai,
217-
workev::Array{$T}, bmat, n, evtype, nev, TOL::Array{$T}, resid::Array{$T}, ncv, v, ldv,
218-
iparam, ipntr, workd::Array{$T}, workl::Array{$T}, lworkl, info)
217+
workev::Vector{$T}, bmat, n, evtype, nev, TOL::Ref{$T}, resid::Vector{$T}, ncv, v, ldv,
218+
iparam, ipntr, workd::Vector{$T}, workl::Vector{$T}, lworkl, info)
219219
ccall(($(string(neupd_name)), :libarpack), Void,
220-
(Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{$T},
221-
Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{$T}, Ptr{UInt8}, Ptr{BlasInt},
222-
Ptr{UInt8}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{$T},
223-
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T},
224-
Ptr{BlasInt}, Ptr{BlasInt}, Clong, Clong, Clong),
225-
&rvec, howmny, select, dr, di, z, &ldz, sigmar, sigmai,
226-
workev, bmat, &n, evtype, &nev, TOL, resid, &ncv, v, &ldv,
227-
iparam, ipntr, workd, workl, &lworkl, info,
228-
sizeof(howmny), sizeof(bmat), sizeof(evtype))
220+
(Ref{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{$T}, Ref{BlasInt},
221+
Ref{$T}, Ref{$T}, Ptr{$T}, Ptr{UInt8}, Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt},
222+
Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ptr{$T}, Ref{BlasInt},
223+
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ref{BlasInt}, Clong, Clong, Clong),
224+
rvec, howmny, select, dr, di, z, ldz,
225+
sigmar, sigmai, workev, bmat, n, evtype, nev,
226+
TOL, resid, ncv, v, ldv,
227+
iparam, ipntr, workd, workl, lworkl, info, sizeof(howmny), sizeof(bmat), sizeof(evtype))
229228
end
230229

231-
function saupd(ido, bmat, n, which, nev, TOL::Array{$T}, resid::Array{$T}, ncv, v::Array{$T}, ldv,
232-
iparam, ipntr, workd::Array{$T}, workl::Array{$T}, lworkl, info)
230+
function saupd(ido, bmat, n, which, nev, TOL::Ref{$T}, resid::Vector{$T}, ncv, v::Matrix{$T}, ldv,
231+
iparam, ipntr, workd::Vector{$T}, workl::Vector{$T}, lworkl, info)
233232
ccall(($(string(saupd_name)), :libarpack), Void,
234-
(Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt},
235-
Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{$T}, Ptr{BlasInt},
236-
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{BlasInt}, Clong, Clong),
237-
ido, bmat, &n, which, &nev, TOL, resid, &ncv, v, &ldv,
238-
iparam, ipntr, workd, workl, &lworkl, info, sizeof(bmat), sizeof(which))
233+
(Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt},
234+
Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ptr{$T}, Ref{BlasInt},
235+
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ref{BlasInt}, Clong, Clong),
236+
ido, bmat, n, which, nev,
237+
TOL, resid, ncv, v, ldv,
238+
iparam, ipntr, workd, workl, lworkl, info, sizeof(bmat), sizeof(which))
239239
end
240240

241241
function seupd(rvec, howmny, select, d, z, ldz, sigma,
242-
bmat, n, evtype, nev, TOL::Array{$T}, resid::Array{$T}, ncv, v::Array{$T}, ldv,
243-
iparam, ipntr, workd::Array{$T}, workl::Array{$T}, lworkl, info)
242+
bmat, n, evtype, nev, TOL::Ref{$T}, resid::Vector{$T}, ncv, v::Matrix{$T}, ldv,
243+
iparam, ipntr, workd::Vector{$T}, workl::Vector{$T}, lworkl, info)
244244
ccall(($(string(seupd_name)), :libarpack), Void,
245-
(Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{$T},
246-
Ptr{UInt8}, Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt},
247-
Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{$T}, Ptr{BlasInt}, Ptr{BlasInt},
248-
Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{BlasInt}, Clong, Clong, Clong),
249-
&rvec, howmny, select, d, z, &ldz, sigma,
250-
bmat, &n, evtype, &nev, TOL, resid, &ncv, v, &ldv,
251-
iparam, ipntr, workd, workl, &lworkl, info, sizeof(howmny), sizeof(bmat), sizeof(evtype))
245+
(Ref{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt},
246+
Ptr{$T}, Ptr{UInt8}, Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt},
247+
Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ptr{$T}, Ref{BlasInt},
248+
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ref{BlasInt}, Clong, Clong, Clong),
249+
rvec, howmny, select, d, z, ldz,
250+
sigma, bmat, n, evtype, nev,
251+
TOL, resid, ncv, v, ldv,
252+
iparam, ipntr, workd, workl, lworkl, info, sizeof(howmny), sizeof(bmat), sizeof(evtype))
252253
end
253254
end
254255
end
@@ -257,30 +258,31 @@ for (T, TR, naupd_name, neupd_name) in
257258
((:Complex128, :Float64, :znaupd_, :zneupd_),
258259
(:Complex64, :Float32, :cnaupd_, :cneupd_))
259260
@eval begin
260-
function naupd(ido, bmat, n, evtype, nev, TOL::Array{$TR}, resid::Array{$T}, ncv, v::Array{$T}, ldv,
261-
iparam, ipntr, workd::Array{$T}, workl::Array{$T}, lworkl,
262-
rwork::Array{$TR}, info)
261+
function naupd(ido, bmat, n, evtype, nev, TOL::Ref{$TR}, resid::Vector{$T}, ncv, v::Matrix{$T}, ldv,
262+
iparam, ipntr, workd::Vector{$T}, workl::Vector{$T}, lworkl,
263+
rwork::Vector{$TR}, info)
263264
ccall(($(string(naupd_name)), :libarpack), Void,
264-
(Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt},
265-
Ptr{$TR}, Ptr{$T}, Ptr{BlasInt}, Ptr{$T}, Ptr{BlasInt},
266-
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt},
267-
Ptr{$TR}, Ptr{BlasInt}),
268-
ido, bmat, &n, evtype, &nev, TOL, resid, &ncv, v, &ldv,
269-
iparam, ipntr, workd, workl, &lworkl, rwork, info)
265+
(Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt},
266+
Ptr{$TR}, Ptr{$T}, Ref{BlasInt}, Ptr{$T}, Ref{BlasInt},
267+
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ptr{$TR}, Ref{BlasInt}),
268+
ido, bmat, n, evtype, nev,
269+
TOL, resid, ncv, v, ldv,
270+
iparam, ipntr, workd, workl, lworkl, rwork, info)
270271
end
271272

272-
function neupd(rvec, howmny, select, d, z, ldz, sigma, workev::Array{$T},
273-
bmat, n, evtype, nev, TOL::Array{$TR}, resid::Array{$T}, ncv, v::Array{$T}, ldv,
274-
iparam, ipntr, workd::Array{$T}, workl::Array{$T}, lworkl,
275-
rwork::Array{$TR}, info)
273+
function neupd(rvec, howmny, select, d, z, ldz, sigma, workev::Vector{$T},
274+
bmat, n, evtype, nev, TOL::Ref{$TR}, resid::Vector{$T}, ncv, v::Matrix{$T}, ldv,
275+
iparam, ipntr, workd::Vector{$T}, workl::Vector{$T}, lworkl,
276+
rwork::Vector{$TR}, info)
276277
ccall(($(string(neupd_name)), :libarpack), Void,
277-
(Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt},
278-
Ptr{$T}, Ptr{$T}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{UInt8}, Ptr{BlasInt},
279-
Ptr{$TR}, Ptr{$T}, Ptr{BlasInt}, Ptr{$T}, Ptr{BlasInt}, Ptr{BlasInt},
280-
Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ptr{BlasInt}, Ptr{$TR}, Ptr{BlasInt}),
281-
&rvec, howmny, select, d, z, &ldz, sigma, workev,
282-
bmat, &n, evtype, &nev, TOL, resid, &ncv, v, &ldv,
283-
iparam, ipntr, workd, workl, &lworkl, rwork, info)
278+
(Ref{BlasInt}, Ptr{UInt8}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt},
279+
Ptr{$T}, Ptr{$T}, Ptr{UInt8}, Ref{BlasInt}, Ptr{UInt8}, Ref{BlasInt},
280+
Ptr{$TR}, Ptr{$T}, Ref{BlasInt}, Ptr{$T}, Ref{BlasInt},
281+
Ptr{BlasInt}, Ptr{BlasInt}, Ptr{$T}, Ptr{$T}, Ref{BlasInt}, Ptr{$TR}, Ref{BlasInt}),
282+
rvec, howmny, select, d, z, ldz,
283+
sigma, workev, bmat, n, evtype, nev,
284+
TOL, resid, ncv, v, ldv,
285+
iparam, ipntr, workd, workl, lworkl, rwork, info)
284286
end
285287
end
286288
end

0 commit comments

Comments
 (0)