Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add Support for Matrices with Int32 and UInt32 as Indices Vectors #43

Closed
RoyiAvital opened this issue Sep 20, 2021 · 5 comments · Fixed by #44
Closed

Add Support for Matrices with Int32 and UInt32 as Indices Vectors #43

RoyiAvital opened this issue Sep 20, 2021 · 5 comments · Fixed by #44

Comments

@RoyiAvital
Copy link

RoyiAvital commented Sep 20, 2021

It seems the code fails when the input sparse matrix has indices with types which are not Int64:

using SparseArrays;
using LimitedLDLFactorizations;

mA = sprand(10, 10, 0.1);
vI, vJ, vK = findnz(mA);
mB = sparse(Int32.(vI), Int32.(vJ), vK, mA.m, mA.n);

lldl(mA); #<! Works
lldl(mB); #<! Fails

I get the following error:

julia> lldl(mB)
ERROR: MethodError: no method matching abspermute!(::Vector{Float64}, ::SubArray{Int32, 1, Vector{Int32}, Tuple{UnitRange{Int64}}, true}, ::Int64)

I think the code should support any integer as indices and any float as value. Wasn't that the design goal?

Performance wise, if the input data is 32 Bit having all other computations in 32 but might assist with performance (Memory and calculation throughput).

For the UInt32 case:

mB = sparse(UInt32.(vI), UInt32.(vJ), vK, mA.m, mA.n);
lldl(mB); #<! Fails

The error:

julia> lldl(mB); #<! Fails
ERROR: MethodError: no method matching amd(::SparseMatrixCSC{Float64, UInt32})

It also fails for UInt64.

Remark: I didn't make the matrix symmetric as it is has no significance for the error.

dpo added a commit that referenced this issue Sep 20, 2021
dpo added a commit that referenced this issue Sep 20, 2021
@dpo dpo closed this as completed in #44 Sep 20, 2021
dpo added a commit that referenced this issue Sep 20, 2021
@RoyiAvital
Copy link
Author

OK, It seems version 0.4 the added support works correctly.
Moreover:

using SparseArrays;
using MatrixDepot;
using LimitedLDLFactorizations;
using BenchmarkTools;

mA = matrixdepot("poisson", 100);
mB = convert(SparseMatrixCSC{Float64, Int32}, mA);

@btime lldl($mA);
@btime lldl($mB);

Gives me ~15% speedup.

@RoyiAvital
Copy link
Author

One more remark, the updated code won't work with UInt32 or UInt64 because amd() doesn't support those (On the C code level).

I wish someone implemented amd() in pure Julia.

@dpo
Copy link
Member

dpo commented Sep 21, 2021

@Sacha0 did at some point, but I can't recall the name or status of his repo.

@dpo
Copy link
Member

dpo commented Sep 21, 2021

Actually here it is: https://github.com/Sacha0/ApproxMinimumDegree.jl

@RoyiAvital
Copy link
Author

That's nice. I wish he continued working on it.
But anyhow, until someone has supported AMD with UInt support there is no real way to have it.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

Successfully merging a pull request may close this issue.

2 participants