Skip to content

Commit 1069f79

Browse files
committed
Update ksizeof for workspaces of block-Krylov solvers
1 parent a2b42ee commit 1069f79

File tree

3 files changed

+22
-22
lines changed

3 files changed

+22
-22
lines changed

docs/src/storage.md

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -149,4 +149,4 @@ Base.format_bytes(free_nbytes) # Total free memory in RAM in bytes.
149149

150150
!!! note
151151
- Beyond having faster operations, using low precisions, such as simple precision, allows to store more coefficients in RAM and solve larger linear problems.
152-
- In the file [test_allocations.jl](https://github.com/JuliaSmoothOptimizers/Krylov.jl/blob/main/test/test_allocations.jl), we use the macro `@allocated` to test that we match the expected storage requirement of each method with a tolerance of 2%.
152+
- In the file [test_allocations.jl](https://github.com/JuliaSmoothOptimizers/Krylov.jl/blob/main/test/test_allocations.jl), we use the macro `@allocated` to verify that the storage requirements of each Krylov solver match the expected values, within a tolerance of 2% (and 5% for block Krylov solvers). These tests are performed across the four main precisions: `Float32`, `Float64`, `ComplexF32`, and `ComplexF64`.

src/krylov_show.jl

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -1,9 +1,9 @@
11
import Base.show, Base.sizeof, Base.format_bytes
22

33
function ksizeof(attribute)
4-
if isa(attribute, Vector{<:AbstractVector}) && !isempty(attribute)
5-
# A vector of vectors is a vector of pointers in Julia.
6-
# All vectors inside a vector have the same size in Krylov.jl
4+
if isa(attribute, Vector{<:AbstractArray}) && !isempty(attribute)
5+
# A vector of arrays is a vector of pointers in Julia.
6+
# All arrays inside a vector have the same size in Krylov.jl
77
size_attribute = sizeof(attribute) + length(attribute) * ksizeof(attribute[1])
88
else
99
size_attribute = sizeof(attribute)

test/test_allocations.jl

Lines changed: 18 additions & 18 deletions
Original file line numberDiff line numberDiff line change
@@ -12,7 +12,7 @@
1212
b = Ao * ones(FC, k) # Dimension m
1313
c = Au * ones(FC, n) # Dimension k
1414
B = A * Matrix{FC}(I, m, p) # Dimension m × p
15-
mem = 200
15+
mem = 100
1616

1717
T = real(FC)
1818
shifts = T[1; 2; 3; 4; 5]
@@ -472,7 +472,7 @@
472472
(x, stats) = lslq(Ao, b) # warmup
473473
actual_lslq_bytes = @allocated lslq(Ao, b)
474474
if VERSION < v"1.11.5" || !Sys.isapple()
475-
@test expected_lslq_bytes actual_lslq_bytes 1.025 * expected_lslq_bytes
475+
@test expected_lslq_bytes actual_lslq_bytes 1.02 * expected_lslq_bytes
476476
end
477477

478478
workspace = LslqWorkspace(Ao, b)
@@ -723,21 +723,21 @@
723723
# - mem (2p*p)-matrices: H
724724
# - lwork-vector: buffer
725725
function storage_block_gmres_bytes(mem, n, p)
726-
res = (2*n*p + p*p + 2p*p + mem*p + mem*n*p + mem*p*p + mem*(mem+1)*p*p/2 + mem*2p*p)
726+
res = (2*n*p + p*p + 2p*p + mem*p + mem*n*p + mem*p*p + (mem*(mem+1)÷2)*p*p + mem*2p*p)
727727
return nbits_FC * res
728728
end
729729

730-
expected_block_gmres_bytes = storage_block_gmres_bytes(mem, n, p)
731-
block_gmres(A, B; memory=mem, itmax=mem) # warmup
732-
actual_block_gmres_bytes = @allocated block_gmres(A, B; memory=mem, itmax=mem)
733-
if VERSION < v"1.11.5" || !Sys.isapple()
734-
@test expected_block_gmres_bytes actual_block_gmres_bytes 1.08 * expected_block_gmres_bytes
735-
end
736-
737730
workspace = BlockGmresWorkspace(A, B; memory=mem)
738731
block_gmres!(workspace, A, B) # warmup
739732
inplace_block_gmres_bytes = @allocated block_gmres!(workspace, A, B)
740733
@test inplace_block_gmres_bytes == 0
734+
735+
expected_block_gmres_bytes = storage_block_gmres_bytes(mem, n, p) + sizeof(workspace.buffer)
736+
block_gmres(A, B; memory=mem, itmax=mem) # warmup
737+
actual_block_gmres_bytes = @allocated block_gmres(A, B; memory=mem, itmax=mem)
738+
if VERSION < v"1.11.5" || !Sys.isapple()
739+
@test expected_block_gmres_bytes actual_block_gmres_bytes 1.05 * expected_block_gmres_bytes
740+
end
741741
end
742742

743743
@testset "BLOCK-MINRES" begin
@@ -752,21 +752,21 @@
752752
# - mem (2p*p)-matrices: H
753753
# - lwork-vector: buffer
754754
function storage_block_minres_bytes(mem, n, p)
755-
res = (2*n*p + p*p + 2p*p + mem*p + mem*n*p + mem*p*p + mem*(mem+1)*p*p/2 + mem*2p*p)
755+
res = (2*n*p + p*p + 2p*p + mem*p + mem*n*p + mem*p*p + (mem*(mem+1)÷2)*p*p + mem*2p*p)
756756
return nbits_FC * res
757757
end
758758

759-
expected_block_minres_bytes = storage_block_minres_bytes(mem, n, p)
759+
workspace = BlockMinresWorkspace(A, B)
760+
block_minres!(workspace, A, B) # warmup
761+
inplace_block_minres_bytes = @allocated block_minres!(workspace, A, B)
762+
@test inplace_block_minres_bytes == 0
763+
764+
expected_block_minres_bytes = storage_block_minres_bytes(mem, n, p) + sizeof(workspace.buffer)
760765
block_minres(A, B) # warmup
761766
# actual_block_minres_bytes = @allocated block_minres(A, B)
762767
# if VERSION < v"1.11.5" || !Sys.isapple()
763-
# @test expected_block_minres_bytes ≤ actual_block_minres_bytes ≤ 1.08 * expected_block_minres_bytes
768+
# @test expected_block_minres_bytes ≤ actual_block_minres_bytes ≤ 1.05 * expected_block_minres_bytes
764769
# end
765-
766-
Workspace = BlockMinresWorkspace(A, B)
767-
block_minres!(Workspace, A, B) # warmup
768-
inplace_block_minres_bytes = @allocated block_minres!(Workspace, A, B)
769-
@test inplace_block_minres_bytes == 0
770770
end
771771
end
772772
end

0 commit comments

Comments
 (0)