From 8ecc2935fc22dbc42f16a0b7827c7064ea73bddc Mon Sep 17 00:00:00 2001 From: Alexis Montoison <35051714+amontoison@users.noreply.github.com> Date: Mon, 20 May 2024 16:20:36 -0400 Subject: [PATCH] Add tests with LinearOperators.jl (#863) * Add tests with LinearOperators.jl * Update block_gmres for complex numbers * Update block_gmres.jl --- .buildkite/pipeline.yml | 14 ++++++++++++++ src/block_gmres.jl | 5 +++-- test/cpu/linear_operators.jl | 16 ++++++++++++++++ 3 files changed, 33 insertions(+), 2 deletions(-) create mode 100644 test/cpu/linear_operators.jl diff --git a/.buildkite/pipeline.yml b/.buildkite/pipeline.yml index 120f1f736..358eb9cac 100644 --- a/.buildkite/pipeline.yml +++ b/.buildkite/pipeline.yml @@ -98,3 +98,17 @@ steps: Pkg.instantiate() include("test/cpu/component_arrays.jl")' timeout_in_minutes: 30 + + - label: "CPUs -- LinearOperators.jl" + plugins: + - JuliaCI/julia#v1: + version: "1.10" + agents: + queue: "juliaecosystem" + command: | + julia --color=yes --project -e ' + using Pkg + Pkg.add("LinearOperators") + Pkg.instantiate() + include("test/cpu/linear_operators.jl")' + timeout_in_minutes: 30 diff --git a/src/block_gmres.jl b/src/block_gmres.jl index c1e04164d..ea0b4b989 100644 --- a/src/block_gmres.jl +++ b/src/block_gmres.jl @@ -148,6 +148,7 @@ kwargs_block_gmres = (:M, :N, :ldiv, :restart, :reorthogonalization, :atol, :rto # Define the blocks D1 and D2 D1 = view(D, 1:p, :) D2 = view(D, p+1:2p, :) + trans = FC <: AbstractFloat ? 'T' : 'C' # Coefficients for mul! α = -one(FC) @@ -263,7 +264,7 @@ kwargs_block_gmres = (:M, :N, :ldiv, :restart, :reorthogonalization, :atol, :rto for i = 1 : inner_iter-1 D1 .= R[nr+i] D2 .= R[nr+i+1] - @kormqr!('L', 'T', H[i], τ[i], D) + @kormqr!('L', trans, H[i], τ[i], D) R[nr+i] .= D1 R[nr+i+1] .= D2 end @@ -276,7 +277,7 @@ kwargs_block_gmres = (:M, :N, :ldiv, :restart, :reorthogonalization, :atol, :rto # Update Zₖ = (Qₖ)ᴴΓE₁ = (Λ₁, ..., Λₖ, Λbarₖ₊₁) D1 .= Z[inner_iter] D2 .= zero(FC) - @kormqr!('L', 'T', H[inner_iter], τ[inner_iter], D) + @kormqr!('L', trans, H[inner_iter], τ[inner_iter], D) Z[inner_iter] .= D1 # Update residual norm estimate. diff --git a/test/cpu/linear_operators.jl b/test/cpu/linear_operators.jl new file mode 100644 index 000000000..3a58c5d20 --- /dev/null +++ b/test/cpu/linear_operators.jl @@ -0,0 +1,16 @@ +using LinearAlgebra, SparseArrays, Test +using Krylov, LinearOperators + +@testset "LinearOperators" begin + n = 50 + p = 5 + + for T in (Float64, ComplexF64) + A = rand(T, n, n) + B = rand(T, n, p) + + opA = LinearOperator(A) + x, stats = block_gmres(opA, B) + @test stats.solved + end +end