diff --git a/Project.toml b/Project.toml index 58a28aba..cce66e6f 100644 --- a/Project.toml +++ b/Project.toml @@ -9,6 +9,7 @@ FastGaussQuadrature = "442a2c76-b920-505d-bb47-c5924d526838" ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210" LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e" LinearSolve = "7ed4a6bd-45f5-4d41-b270-4a48e9bafcae" +NearestNeighbors = "b8a86587-4115-5ab1-83bc-aa920d37bbce" Printf = "de0858da-6303-5e67-8744-51eddeeeb8d7" SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf" StaticArrays = "90137ffa-7385-5640-81b9-e52037218182" diff --git a/src/Bcube.jl b/src/Bcube.jl index 58c9763c..72d46cfd 100644 --- a/src/Bcube.jl +++ b/src/Bcube.jl @@ -13,6 +13,7 @@ using Printf # just for tmp vtk, to be removed # import LinearSolve: solve, solve!, LinearProblem import LinearSolve using Symbolics # used for generation of Lagrange shape functions +using NearestNeighbors const MAX_LENGTH_STATICARRAY = (10^6) @@ -96,6 +97,8 @@ include("./mapping/mapping.jl") include("./mapping/ref2phys.jl") export get_cell_centers +include("mapping/findpoint.jl") + include("./cellfunction/eval_point.jl") include("./cellfunction/cellfunction.jl") diff --git a/src/mapping/findpoint.jl b/src/mapping/findpoint.jl new file mode 100644 index 00000000..401ab20a --- /dev/null +++ b/src/mapping/findpoint.jl @@ -0,0 +1,47 @@ +struct PointFinder{T, M} + tree::T + mesh::M +end + +function PointFinder(mesh::AbstractMesh) + xc0 = get_cell_centers(mesh) + xc = reshape( + [xc0[n][idim] for n in 1:length(xc0) for idim in 1:spacedim(mesh)], + spacedim(mesh), + length(xc0), + ) + tree = KDTree(xc; leafsize = 10) + PointFinder{typeof(tree), typeof(mesh)}(tree, mesh) +end + +function find_cell(pf::PointFinder, point) + idxs, dists = knn(pf.tree, point, 1) + c2c_n = connectivity_cell2cell_by_nodes(pf.mesh) + cellids = [idxs[1], c2c_n[idxs[1]]...] + + for icell in cellids + cinfo = CellInfo(pf.mesh, icell) + cpoint = CellPoint(point, cinfo, PhysicalDomain()) + isIn = point_in_cell(cinfo, cpoint) + isIn && return icell + end + return nothing +end + +function point_in_cell(cinfo, cpoint) + cpoint_ref = change_domain(cpoint, ReferenceDomain()) + point_in_shape(shape(celltype(cinfo)), get_coords(cpoint_ref)) +end + +function point_in_shape(s::Square, x) + get_coords(s)[1][1] ≤ x[1] ≤ get_coords(s)[3][1] && + get_coords(s)[1][2] ≤ x[2] ≤ get_coords(s)[3][2] +end + +function point_in_shape(shape::AbstractShape, x) + for (normal, f2n) in zip(normals(shape), faces2nodes(shape)) + dx = (x - get_coords(shape)[first(f2n)]) + (dx ⋅ normal > 0) && (return false) + end + return true +end \ No newline at end of file diff --git a/test/interpolation/test_cellfunction.jl b/test/interpolation/test_cellfunction.jl index d0783602..a3312b0a 100644 --- a/test/interpolation/test_cellfunction.jl +++ b/test/interpolation/test_cellfunction.jl @@ -287,4 +287,34 @@ end @test v1 ⋅ u ≈ v1_in_2 ⋅ u @test abs(ν2 ⋅ v1_in_2) < 1e-16 end + + @testset "Point interpolation" begin + degree = 2 + path = joinpath(tempdir, "mesh.msh") + gen_rectangle_mesh_with_tri_and_quad(path; nx = 10, ny = 10, xc = 0.5, yc = 0.5) + mesh = read_msh(path) + Uspace = TrialFESpace(FunctionSpace(:Lagrange, degree), mesh; size = 2) + Vspace = TestFESpace(Uspace) + dΩ = Measure(CellDomain(mesh), 2 * degree + 1) + + pointFinder = Bcube.PointFinder(mesh) + + u = FEFunction(Uspace) + f1((x, y)) = SA[x + 4y + 1, 7x * x + 2y * x + 2y * y + 3y + 12] + projection_l2!(u, PhysicalFunction(f1), dΩ) + npoints = 20 + xp = [rand(SVector{2}) for i in 1:npoints] + icells = map(Base.Fix1(Bcube.find_cell, pointFinder), xp) + cpoints = map( + (x, i) -> + Bcube.CellPoint(x, Bcube.CellInfo(mesh, i), Bcube.PhysicalDomain()), + xp, + icells, + ) + up = map(cpoints) do cpoint + uᵢ = Bcube.materialize(u, Bcube.get_cellinfo(cpoint)) + Bcube.materialize(uᵢ, cpoint) + end + @test all(f1.(xp) .≈ up) + end end