Skip to content

Commit

Permalink
add tests for blocked solve
Browse files Browse the repository at this point in the history
  • Loading branch information
PetrKryslUCSD committed Oct 12, 2023
1 parent 7614e83 commit 0bdf551
Show file tree
Hide file tree
Showing 3 changed files with 95 additions and 10 deletions.
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "FinEtools"
uuid = "91bb5406-6c9a-523d-811d-0644c4229550"
authors = ["Petr Krysl <[email protected]>"]
version = "7.1.3"
version = "7.1.4"

[deps]
ChunkSplitters = "ae650224-84b6-46f8-82ea-d812ca08434e"
Expand Down
28 changes: 21 additions & 7 deletions src/AlgoBaseModule.jl
Original file line number Diff line number Diff line change
Expand Up @@ -629,16 +629,30 @@ function vector_blocked(V, nfreedofs)
end

"""
solve!(u::F, K::M, F::V) where {F<:AbstractField, M<:AbstractMatrix, V<:AbstractVector}
solve_blocked(A::M, b::VB, x::VX, nfreedofs::IT) where {M<:AbstractMatrix, VB<:AbstractVector, VX<:AbstractVector, IT<:Integer}
Solve a blocked system of linear algebraic equations.
b_f and x_d are known, x_f and b_d need to be computed.
"""
function solve_blocked(A::M, b::VB, x::VX, nfreedofs::IT) where {M<:AbstractMatrix, VB<:AbstractVector, VX<:AbstractVector, IT<:Integer}
A_b = matrix_blocked(A, nfreedofs)
b_b = vector_blocked(b, nfreedofs)
x_b = vector_blocked(x, nfreedofs)
x_f = A_b.ff \ (b_b.f - A_b.fd * x_b.d)
b_d = A_b.df * x_f + A_b.dd * x_b.d
return x_f, b_d
end

"""
solve_blocked!(u::AF, K::M, F::V) where {AF<:AbstractField, M<:AbstractMatrix, V<:AbstractVector}
Solve a system of linear algebraic equations.
"""
function solve!(u::U, K::M, F::V) where {U<:AbstractField, M<:AbstractMatrix, V<:AbstractVector}
K_ff, K_fd = matrix_blocked(K, nfreedofs(u), nfreedofs(u), (:ff, :fd))[(:ff, :fd)]
F_f = vector_blocked(F, nfreedofs(u))[:f]
U_d = gathersysvec(u, :d)
U_f = K_ff\(F_f - K_fd * U_d)
scattersysvec!(u, U_f)
function solve_blocked!(u::AF, K::M, F::V) where {AF<:AbstractField, M<:AbstractMatrix, V<:AbstractVector}
U = gathersysvec(u, :a)
x_f, b_d = solve_blocked(K, F, U, nfreedofs(u))
scattersysvec!(u, x_f)
return u
end

Expand Down
75 changes: 73 additions & 2 deletions test/test_basics.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2527,7 +2527,7 @@ end

module mblocked1
using Test
using FinEtools.AlgoBaseModule: solve!, matrix_blocked, vector_blocked
using FinEtools.AlgoBaseModule: matrix_blocked, vector_blocked
using LinearAlgebra

function test(n, nfreedofs)
Expand Down Expand Up @@ -2575,7 +2575,7 @@ end

module mblocked2
using Test
using FinEtools.AlgoBaseModule: solve!, vector_blocked
using FinEtools.AlgoBaseModule: vector_blocked
using LinearAlgebra

function test(n, nfreedofs)
Expand Down Expand Up @@ -2617,3 +2617,74 @@ nothing
end


module mblocked3
using Test
using FinEtools.AlgoBaseModule: solve_blocked, solve_blocked!, vector_blocked, matrix_blocked
using LinearAlgebra

function test(n, nfreedofs)
A = rand(n, n) + I(n)
b = rand(n)
x = A \ b

b_b = vector_blocked(b, nfreedofs)
x_b = vector_blocked(x, nfreedofs)
A_b = matrix_blocked(A, nfreedofs)

x_f = A_b.ff \ (b_b.f - A_b.fd * x_b.d)
b_d = A_b.df * x_f + A_b.dd * x_b.d

@test norm(x - vcat(x_f, x_b.d)) / norm(x) < 1.0e-9
@test norm(x_f - x_b.f) / norm(x_f) < 1.0e-9
@test norm(b_b.d - b_d) / norm(b_d) < 1.0e-9

x_f, b_d = solve_blocked(A, b, x, nfreedofs)

@test norm(x - vcat(x_f, x_b.d)) / norm(x) < 1.0e-9
@test norm(x_f - x_b.f) / norm(x_f) < 1.0e-9
@test norm(b_b.d - b_d) / norm(b_d) < 1.0e-9

true
end
test(10, 9)
test(100, 9)
test(1000, 999)
test(200, 99)
nothing
end

module mblocked4
using Test
using FinEtools
using FinEtools.AlgoBaseModule: solve_blocked, solve_blocked!, vector_blocked, matrix_blocked
using LinearAlgebra

function test(n, nfreed)
A = rand(n, n) + 2*I(n)
b = rand(n)
x = A \ b

u = NodalField(x)
u.is_fixed[nfreed+1:end] .= true
numberdofs!(u)

solve_blocked!(u, A, b)

x_f = gathersysvec(u, :f)
x_b = vector_blocked(x, nfreed)

@test norm(x_f - x_b.f) / norm(x_f) < 1.0e-9

true
end
test(10, 9)
test(100, 9)
test(1000, 999)
test(200, 99)
test(130, 99)
test(130, 1)
nothing
end



2 comments on commit 0bdf551

@PetrKryslUCSD
Copy link
Owner Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@JuliaRegistrator
Copy link

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Registration pull request created: JuliaRegistries/General/93268

After the above pull request is merged, it is recommended that a tag is created on this repository for the registered package version.

This will be done automatically if the Julia TagBot GitHub Action is installed, or can be done manually through the github interface, or via:

git tag -a v7.1.4 -m "<description of version>" 0bdf55137ac29a37a6d3282a4d999b86f4a237e9
git push origin v7.1.4

Please sign in to comment.