From 62a929fbdab30a8afa080419328636f5288f6abb Mon Sep 17 00:00:00 2001 From: mbouyges Date: Thu, 28 Mar 2024 17:13:23 +0100 Subject: [PATCH] compute returns SparseVector --- src/assembler/assembler.jl | 12 ++++++++---- test/Project.toml | 3 ++- test/dof/test_assembler.jl | 8 ++++---- test/integration/test_integration.jl | 21 +++++++++++---------- test/operator/test_algebra.jl | 3 +-- test/runtests.jl | 1 + 6 files changed, 27 insertions(+), 21 deletions(-) diff --git a/src/assembler/assembler.jl b/src/assembler/assembler.jl index 001f3ca2..9c9bc801 100644 --- a/src/assembler/assembler.jl +++ b/src/assembler/assembler.jl @@ -734,17 +734,18 @@ Base.repeat(a::SVector, ::Val{N}) where {N} = reduce(vcat, ntuple(i -> a, Val(N) Compute an integral, independently from a FEM/DG framework (i.e without FESpace) -Return a dict (same size of the domain of integration) => . +Return a `SparseVector`. The indices of the domain elements are used to store +the result of the integration in this sparse vector. # Example Compute volume of each cell and each face. ```julia mesh = rectangle_mesh(2, 3) dΩ = Measure(CellDomain(mesh), 1) -∂Ω = Measure(BoundaryFaceDomain(mesh), 1) +dΓ = Measure(BoundaryFaceDomain(mesh), 1) f = PhysicalFunction(x -> 1) @show Bcube.compute(∫(f)dΩ) -@show Bcube.compute(∫(side⁻(f))∂Ω) +@show Bcube.compute(∫(side⁻(f))dΓ) ``` """ function compute(integration::Integration) @@ -758,9 +759,12 @@ function compute(integration::Integration) _integrate_on_ref_element(_f, elementInfo, quadrature) end - return Dict(indices(domain) .=> values) + return SparseVector(_domain_to_mesh_nelts(domain), collect(indices(domain)), values) end +_domain_to_mesh_nelts(domain::AbstractCellDomain) = ncells(get_mesh(domain)) +_domain_to_mesh_nelts(domain::AbstractFaceDomain) = nfaces(get_mesh(domain)) + """ AbstractFaceSidePair{A} <: AbstractLazyWrap{A} diff --git a/test/Project.toml b/test/Project.toml index b176adf1..56da7382 100644 --- a/test/Project.toml +++ b/test/Project.toml @@ -5,4 +5,5 @@ StaticArrays = "90137ffa-7385-5640-81b9-e52037218182" LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e" DelimitedFiles = "8bb1440f-4735-579b-a4ab-409b98df4dab" SHA = "ea8e919c-243c-51af-8825-aaa63cd721ce" -ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210" \ No newline at end of file +ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210" +SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf" \ No newline at end of file diff --git a/test/dof/test_assembler.jl b/test/dof/test_assembler.jl index 78cea2d8..6502c870 100644 --- a/test/dof/test_assembler.jl +++ b/test/dof/test_assembler.jl @@ -264,8 +264,8 @@ sys = Bcube.AffineFESystem(a, l, U, V) uh = Bcube.solve(sys) - l2(u) = sqrt(sum(values(compute(∫(u ⋅ u)dΩ)))) - h1(u) = sqrt(sum(values(compute(∫(u ⋅ u + ∇(u) ⋅ ∇(u))dΩ)))) + l2(u) = sqrt(sum(compute(∫(u ⋅ u)dΩ))) + h1(u) = sqrt(sum(compute(∫(u ⋅ u + ∇(u) ⋅ ∇(u))dΩ))) e = uref - uh el2 = l2(e) @@ -313,7 +313,7 @@ int_val = 4.0 / 3.0 - volume = sum(values(compute(∫(PhysicalFunction(x -> 1.0))dΩ))) + volume = sum(compute(∫(PhysicalFunction(x -> 1.0))dΩ)) # Define bilinear and linear forms function a((u, λᵤ), (v, λᵥ)) @@ -339,7 +339,7 @@ u_ref = PhysicalFunction(x -> -0.5 * (x[1] - 1.0)^2 + 1.0) error = u_ref - u - l2(u) = sqrt(sum(values(compute(∫(u ⋅ u)dΩ)))) + l2(u) = sqrt(sum(compute(∫(u ⋅ u)dΩ))) el2 = l2(error) tol = 1.e-15 @test el2 < tol diff --git a/test/integration/test_integration.jl b/test/integration/test_integration.jl index 9128c8ec..bff1d33c 100644 --- a/test/integration/test_integration.jl +++ b/test/integration/test_integration.jl @@ -381,7 +381,7 @@ end Γ = BoundaryFaceDomain(mesh, "boundary") dΓ = Measure(Γ, 2) nΓ = get_face_normals(dΓ) - val = sum(values(compute(∫(side_n(u) ⋅ side_n(nΓ))dΓ))) + val = sum(compute(∫(side_n(u) ⋅ side_n(nΓ))dΓ)) @test isapprox(val, 0.0, atol = 1e-15) end @@ -420,7 +420,7 @@ end Γ = BoundaryFaceDomain(mesh, "boundary") dΓ = Measure(Γ, 2) nΓ = get_face_normals(dΓ) - val = sum(values(compute(∫(side_n(u) ⋅ side_n(nΓ))dΓ))) + val = sum(compute(∫(side_n(u) ⋅ side_n(nΓ))dΓ)) @test isapprox(val, 0.0, atol = 1e-14) end @@ -463,7 +463,7 @@ end Γ = BoundaryFaceDomain(mesh, "boundary") dΓ = Measure(Γ, 2) nΓ = get_face_normals(dΓ) - val = sum(values(compute(∫(side_n(u) ⋅ side_n(nΓ))dΓ))) + val = sum(compute(∫(side_n(u) ⋅ side_n(nΓ))dΓ)) @test isapprox(val, 0.0, atol = 1e-15) end @@ -471,7 +471,7 @@ end mesh = translate(one_cell_mesh(:line), [1.0]) dΩ = Measure(CellDomain(mesh), 2) g = PhysicalFunction(x -> 2 * x[1]) - b = collect(values(compute(∫(g)dΩ))) + b = compute(∫(g)dΩ) @test b[1] == 4.0 end @@ -521,10 +521,10 @@ end @test Bcube.integrate(x -> 1, ctype, cnodes, Quadrature(1)) ≈ 4 dΩ = Measure(CellDomain(mesh), 2) g = PhysicalFunction(x -> 1) - b = collect(values(compute(∫(g)dΩ))) + b = compute(∫(g)dΩ) @test b[1] ≈ 4 gref = ReferenceFunction(x -> 1) - b = collect(values(compute(∫(gref)dΩ))) + b = compute(∫(gref)dΩ) @test b[1] ≈ 4 lx = 2.0 @@ -548,8 +548,9 @@ end surface_ref = (s(-e12, e23), s(-e12, e24), s(e24, e23), s(e13, e14)) for i in 1:4 dΓ = Measure(BoundaryFaceDomain(mesh, "F$i"), 1) - b = collect(values(compute(∫(side⁻(g))dΓ))) - @test b[1] ≈ surface_ref[i] + b = compute(∫(side⁻(g))dΓ) + I, = findnz(b) + @test b[I[1]] ≈ surface_ref[i] end end @@ -571,7 +572,7 @@ end @test Bcube.integrate(x -> 1, ctype, cnodes, Quadrature(1)) == 0.75 dΩ = Measure(CellDomain(mesh), 2) g = PhysicalFunction(x -> 1) - b = collect(values(compute(∫(g)dΩ))) + b = compute(∫(g)dΩ) @test b[1] ≈ 0.75 # Whole cylinder : build a cylinder of radius 1 and length 1, and compute its volume @@ -581,7 +582,7 @@ end dΩ = Measure(CellDomain(mesh), 2) g = PhysicalFunction(x -> 1) - b = values(compute(∫(g)dΩ)) + b = compute(∫(g)dΩ) @test sum(b) ≈ 3.1365484905459 end end diff --git a/test/operator/test_algebra.jl b/test/operator/test_algebra.jl index 0f4a1c86..79e0b292 100644 --- a/test/operator/test_algebra.jl +++ b/test/operator/test_algebra.jl @@ -55,8 +55,7 @@ translate!(mesh, SA[-0.5, 1.0]) # the translation vector can be anything scale!(mesh, 2.0) dΩ = Measure(CellDomain(mesh), 1) - @test collect(values(compute(∫(∇(PhysicalFunction(x -> x[1])) ⋅ [1, 1])dΩ)))[1] == - 16.0 + @test compute(∫(∇(PhysicalFunction(x -> x[1])) ⋅ [1, 1])dΩ)[1] == 16.0 # Gradient of a vector PhysicalFunction mesh = one_cell_mesh(:quad) diff --git a/test/runtests.jl b/test/runtests.jl index 020f930d..d68160ed 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -5,6 +5,7 @@ using LinearAlgebra using DelimitedFiles using ForwardDiff using SHA: sha1 +using SparseArrays # from : # https://discourse.julialang.org/t/what-general-purpose-commands-do-you-usually-end-up-adding-to-your-projects/4889