diff --git a/.github/workflows/test-app.yml b/.github/workflows/test-app.yml index edae7a26..986fc4fa 100644 --- a/.github/workflows/test-app.yml +++ b/.github/workflows/test-app.yml @@ -46,9 +46,12 @@ jobs: steps: - uses: actions/checkout@v2 + # Workaround until https://github.com/h5py/h5py/pull/2101 is merged - name: Install meshio/h5py run: | - CC=mpicc pip3 install --no-cache-dir --no-binary=h5py h5py meshio + pip3 install mpi4py --no-cache-dir --upgrade + HDF5_MPI="ON" HDF5_DIR="/usr/local" pip3 install --no-cache-dir --no-binary=h5py git+https://github.com/julesghub/h5py.git@py310-compat-2100 --upgrade + pip3 install meshio --no-cache-dir --upgrade - name: Get Basix and install uses: actions/checkout@v2 diff --git a/cpp/CMakeLists.txt b/cpp/CMakeLists.txt index 9b1c182d..92d15f5d 100644 --- a/cpp/CMakeLists.txt +++ b/cpp/CMakeLists.txt @@ -55,7 +55,7 @@ target_link_libraries(dolfinx_contact PUBLIC Basix::basix) target_link_libraries(dolfinx_contact PUBLIC dolfinx) include(GNUInstallDirs) -install(FILES Contact.h contact_kernels.hpp utils.h coefficients.h geometric_quantities.h SubMesh.h QuadratureRule.h DESTINATION ${CMAKE_INSTALL_INCLUDEDIR}/dolfinx_contact COMPONENT Development) +install(FILES Contact.h contact_kernels.hpp utils.h coefficients.h geometric_quantities.h SubMesh.h QuadratureRule.h RayTracing.h DESTINATION ${CMAKE_INSTALL_INCLUDEDIR}/dolfinx_contact COMPONENT Development) target_sources(dolfinx_contact PRIVATE ${CMAKE_CURRENT_SOURCE_DIR}/utils.cpp @@ -64,6 +64,7 @@ target_sources(dolfinx_contact PRIVATE ${CMAKE_CURRENT_SOURCE_DIR}/SubMesh.cpp ${CMAKE_CURRENT_SOURCE_DIR}/QuadratureRule.cpp ${CMAKE_CURRENT_SOURCE_DIR}/Contact.cpp + ${CMAKE_CURRENT_SOURCE_DIR}/RayTracing.cpp ) # Set target include location (for build and installed) diff --git a/cpp/RayTracing.cpp b/cpp/RayTracing.cpp new file mode 100644 index 00000000..f4dedd8b --- /dev/null +++ b/cpp/RayTracing.cpp @@ -0,0 +1,288 @@ + +// Copyright (C) 2022 Jørgen S. Dokken +// +// This file is part of DOLFINx_CONTACT +// +// SPDX-License-Identifier: MIT + +#include "RayTracing.h" +#include +#include +#include +#include + +namespace +{ + +/// Get function that parameterizes a facet of a given cell +/// @param[in] cell_type The cell type +/// @param[in] facet_index The facet index (local to cell) +/// @returns Function that computes the coordinate parameterization of the local +/// facet on the reference cell. +std::function>( + xt::xtensor_fixed>)> +get_3D_parameterization(dolfinx::mesh::CellType cell_type, int facet_index) +{ + switch (cell_type) + { + case dolfinx::mesh::CellType::tetrahedron: + break; + case dolfinx::mesh::CellType::hexahedron: + break; + default: + throw std::invalid_argument("Unsupported cell type"); + break; + } + + const int tdim = dolfinx::mesh::cell_dim(cell_type); + + const int num_facets = dolfinx::mesh::cell_num_entities(cell_type, 2); + if (facet_index >= num_facets) + throw std::invalid_argument( + "Invalid facet index (larger than number of facets"); + + // Get basix geometry information + basix::cell::type basix_cell + = dolfinx::mesh::cell_type_to_basix_type(cell_type); + const xt::xtensor x = basix::cell::geometry(basix_cell); + const std::vector> facets + = basix::cell::topology(basix_cell)[tdim - 1]; + + // Create parameterization function exploiting that the mapping between + // reference geometries are affine + std::function>( + xt::xtensor_fixed>)> + func = [x, facet = facets[facet_index]]( + xt::xtensor_fixed> xi) + -> xt::xtensor_fixed> + { + auto x0 = xt::row(x, facet[0]); + xt::xtensor_fixed> vals = x0; + + for (std::size_t i = 0; i < 3; ++i) + for (std::size_t j = 0; j < 2; ++j) + vals(0, i) += (xt::row(x, facet[j + 1])[i] - x0[i]) * xi[j]; + return vals; + }; + return func; +} +//------------------------------------------------------------------------------------------------ +/// Get derivative of the parameterization with respect to the input +/// parameters +/// @param[in] cell_type The cell type +/// @param[in] facet_index The facet index (local to cell) +/// @returns The Jacobian of the parameterization +xt::xtensor_fixed> +get_parameterization_jacobian(dolfinx::mesh::CellType cell_type, + int facet_index) +{ + switch (cell_type) + { + case dolfinx::mesh::CellType::tetrahedron: + break; + case dolfinx::mesh::CellType::hexahedron: + break; + default: + throw std::invalid_argument("Unsupported cell type"); + break; + } + + basix::cell::type basix_cell + = dolfinx::mesh::cell_type_to_basix_type(cell_type); + xt::xtensor facet_jacobians + = basix::cell::facet_jacobians(basix_cell); + + xt::xtensor_fixed> output; + output = xt::view(facet_jacobians, facet_index, xt::all(), xt::all()); + return output; +} + +} // namespace +//------------------------------------------------------------------------------------------------ +int dolfinx_contact::allocated_3D_ray_tracing( + dolfinx_contact::newton_3D_storage& storage, + xt::xtensor& basis_values, xt::xtensor& dphi, + int max_iter, double tol, const dolfinx::fem::CoordinateElement& cmap, + dolfinx::mesh::CellType cell_type, + const xt::xtensor& coordinate_dofs, + const std::function>( + xt::xtensor_fixed>)>& reference_map) +{ + int status = -1; + constexpr int tdim = 3; + storage.x_k = {{0, 0, 0}}; + storage.xi_k = {0.5, 0.25}; + for (int k = 0; k < max_iter; ++k) + { + // Evaluate reference coordinate at current iteration + storage.X_k = reference_map(storage.xi_k); + + // Tabulate coordinate element basis function + cmap.tabulate(1, storage.X_k, basis_values); + + // Push forward reference coordinate + cmap.push_forward(storage.x_k, coordinate_dofs, + xt::view(basis_values, 0, xt::all(), xt::all(), 0)); + dphi = xt::view(basis_values, xt::xrange(1, tdim + 1), 0, xt::all(), 0); + + // Compute Jacobian + std::fill(storage.J.begin(), storage.J.end(), 0); + cmap.compute_jacobian(dphi, coordinate_dofs, storage.J); + + // Compute residual at current iteration + std::fill(storage.Gk.begin(), storage.Gk.end(), 0); + for (std::size_t i = 0; i < 3; ++i) + { + storage.Gk[0] + += (storage.x_k(0, i) - storage.point[i]) * storage.tangents(0, i); + storage.Gk[1] + += (storage.x_k(0, i) - storage.point[i]) * storage.tangents(1, i); + } + + // Check for convergence in first iteration + if ((k == 0) and (std::abs(storage.Gk[0]) < tol) + and (std::abs(storage.Gk[1]) < tol)) + break; + + /// Compute dGk/dxi + std::fill(storage.dGk_tmp.begin(), storage.dGk_tmp.end(), 0); + dolfinx::math::dot(storage.J, storage.dxi, storage.dGk_tmp); + std::fill(storage.dGk.begin(), storage.dGk.end(), 0); + + for (std::size_t i = 0; i < 2; ++i) + for (std::size_t j = 0; j < 2; ++j) + for (std::size_t l = 0; l < 3; ++l) + storage.dGk(i, j) += storage.dGk_tmp(l, j) * storage.tangents(i, l); + + // Invert dGk/dxi + double det_dGk = dolfinx::math::det(storage.dGk); + if (std::abs(det_dGk) < tol) + { + status = -2; + break; + } + dolfinx::math::inv(storage.dGk, storage.dGk_inv); + + // Compute dxi + std::fill(storage.dxi_k.begin(), storage.dxi_k.end(), 0); + for (std::size_t i = 0; i < 2; ++i) + for (std::size_t j = 0; j < 2; ++j) + storage.dxi_k[i] += storage.dGk_inv(i, j) * storage.Gk[j]; + // Check for convergence + if ((storage.dxi_k[0] * storage.dxi_k[0] + + storage.dxi_k[1] * storage.dxi_k[1]) + < tol * tol) + { + status = 1; + break; + } + + // Update xi + std::transform(storage.xi_k.cbegin(), storage.xi_k.cend(), + storage.dxi_k.cbegin(), storage.xi_k.begin(), + [](auto x, auto y) { return x - y; }); + } + // Check if converged parameters are valid + switch (cell_type) + { + + case dolfinx::mesh::CellType::tetrahedron: + if ((storage.xi_k[0] < -tol) or (storage.xi_k[0] > 1 + tol) + or (storage.xi_k[1] < -tol) + or (storage.xi_k[1] > 1 - storage.xi_k[0] + tol)) + { + status = -3; + } + break; + case dolfinx::mesh::CellType::hexahedron: + if ((storage.xi_k[0] < -tol) or (storage.xi_k[0] > 1 + tol) + or (storage.xi_k[1] < -tol) or (storage.xi_k[1] > 1 + tol)) + { + status = -3; + } + break; + default: + throw std::invalid_argument("Unsupported cell type"); + } + return status; +} +//------------------------------------------------------------------------------------------------ +std::tuple>> +dolfinx_contact::compute_3D_ray( + const dolfinx::mesh::Mesh& mesh, + const xt::xtensor_fixed>& point, + const xt::xtensor_fixed>& tangents, + const std::vector>& cells, const int max_iter, + const double tol) +{ + int status = -1; + dolfinx::mesh::CellType cell_type = mesh.topology().cell_type(); + const int tdim = mesh.topology().dim(); + + const dolfinx::fem::CoordinateElement& cmap = mesh.geometry().cmap(); + + // Get cell coordinates/geometry + const dolfinx::mesh::Geometry& geometry = mesh.geometry(); + const dolfinx::graph::AdjacencyList& x_dofmap + = geometry.dofmap(); + const int gdim = geometry.dim(); + xtl::span x_g = geometry.x(); + const std::size_t num_dofs_g = cmap.dim(); + xt::xtensor coordinate_dofs({num_dofs_g, 3}); + xt::xtensor dphi({(std::size_t)tdim, num_dofs_g}); + + if ((gdim != tdim) or (gdim != 3)) + { + throw std::invalid_argument("This raytracing algorithm is specialized " + "for meshes with topological " + "and geometrical dimension 3"); + } + + // Temporary variables + const std::array basis_shape = cmap.tabulate_shape(1, 1); + xt::xtensor basis_values(basis_shape); + + std::size_t cell_idx = -1; + dolfinx_contact::newton_3D_storage allocated_memory; + allocated_memory.tangents = tangents; + allocated_memory.point = point; + for (std::size_t c = 0; c < cells.size(); ++c) + { + + // Get cell geometry + auto [cell, facet_index] = cells[c]; + auto x_dofs = x_dofmap.links(cell); + for (std::size_t j = 0; j < x_dofs.size(); ++j) + { + dolfinx::common::impl::copy_N<3>( + std::next(x_g.begin(), 3 * x_dofs[j]), + std::next(coordinate_dofs.begin(), 3 * j)); + } + // Assign Jacobian of reference mapping + allocated_memory.dxi + = get_parameterization_jacobian(cell_type, facet_index); + // Get parameterization map + auto reference_map = get_3D_parameterization(cell_type, facet_index); + + status = dolfinx_contact::allocated_3D_ray_tracing( + allocated_memory, basis_values, dphi, max_iter, tol, cmap, cell_type, + coordinate_dofs, reference_map); + if (status > 0) + { + cell_idx = c; + break; + } + } + if (status < 0) + LOG(WARNING) << "No ray through the facets have been found"; + + xt::xtensor_fixed> output_coords; + std::copy(allocated_memory.x_k.cbegin(), allocated_memory.x_k.cend(), + output_coords.begin()); + std::copy(allocated_memory.X_k.cbegin(), allocated_memory.X_k.cend(), + std::next(output_coords.begin(), 3)); + + std::tuple>> + output = std::make_tuple(status, cell_idx, output_coords); + return output; +}; \ No newline at end of file diff --git a/cpp/RayTracing.h b/cpp/RayTracing.h new file mode 100644 index 00000000..cdc2ab40 --- /dev/null +++ b/cpp/RayTracing.h @@ -0,0 +1,100 @@ + +// Copyright (C) 2021 Jørgen S. Dokken +// +// This file is part of DOLFINx_CONTACT +// +// SPDX-License-Identifier: MIT + +#pragma once + +#include +#include + +namespace dolfinx_contact +{ +struct newton_3D_storage +{ + xt::xtensor_fixed> + dxi; // Jacobian of reference mapping + xt::xtensor + X_k; // Solution on reference domain (for Newton solver) + xt::xtensor x_k; // Solution in physical space (for Newton solver) + xt::xtensor_fixed> + xi_k; // Reference parameters (for Newton solver) + xt::xtensor_fixed> + dxi_k; // Gradient of reference parameters (for Newton Solver) + xt::xtensor_fixed> J; // Jacobian of the cell + xt::xtensor_fixed> + dGk_tmp; // Temporary variable to invert Jacobian of Newton solver LHS + xt::xtensor_fixed> dGk; // Newton solver LHS Jacobian + xt::xtensor_fixed> + dGk_inv; // Inverse of Newton solver LHS Jacobian + xt::xtensor_fixed> + Gk; // Residual (RHS) of Newton solver + xt::xtensor_fixed> tangents; // Tangents of ray + xt::xtensor_fixed> point; // Point of origin for ray +}; + +/// @brief Compute raytracing with no dynamic memory allocation +/// +/// Implementation of 3D ray-tracing, using no dynamic memory allocation +/// +/// @param[in,out] storage Structure holding all memory required for +/// the newton iteration. +/// @note It is expected that the variables tangents, point, xi is filled with +/// appropriate values +/// @param[in, out] basis_values Four-dimensional array to write basis values +/// into. +/// @param[in, out] dphi Two-dimensional matrix to write the derviative of the +/// basis functions into +/// @param[in] max_iter Maximum number of iterations for the Newton solver +/// @param[in] tol The tolerance for termination of the Newton solver +/// @param[in] cmap The coordinate element +/// @param[in] cell_type The cell type of the mesh +/// @param[in] coordinate_dofs The cell geometry +/// @param[in] reference_map Function mapping from reference parameters (xi, +/// eta) to the physical element +int allocated_3D_ray_tracing( + newton_3D_storage& storage, xt::xtensor& basis_values, + xt::xtensor& dphi, int max_iter, double tol, + const dolfinx::fem::CoordinateElement& cmap, + dolfinx::mesh::CellType cell_type, + const xt::xtensor& coordinate_dofs, + const std::function>( + xt::xtensor_fixed>)>& reference_map); + +/// @brief Compute the intersection between a ray and a facet in the mesh. +/// +/// The implementation solves dot(\Phi(\xi, \eta)-p, t_i)=0, i=1,2 +/// where \Phi(\xi,\eta) is the parameterized mapping from the reference +/// facet to the physical facet, p the point of origin of the ray, and t_1, +/// t_2 the two tangents defining the ray. For more details, see +/// DOI: 10.1016/j.compstruc.2015.02.027 (eq 14). +/// +/// @note The problem is solved using Newton's method +/// +/// @param[in] mesh The mesh +/// @param[in] point The point of origin for the ray +/// @param[in] tangents The tangents of the ray. Each row corresponds to a +/// tangent. +/// @param[in] cells List of tuples (cell, facet), where the cell index is +/// local to process and the facet index is local to the cell cell +/// @param[in] max_iter The maximum number of iterations to use for Newton's +/// method +/// @param[in] tol The tolerance for convergence in Newton's method +/// @returns A triplet (status, cell_idx, points), where x is the +/// convergence status, cell_idx is which entry in the input list the ray +/// goes through and point is the collision point in global and reference +/// coordinates. +/// @note The convergence status is >0 if converging, -1 if the facet is if +/// the maximum number of iterations are reached, -2 if the facet is +/// parallel with the tangent, -3 if the Newton solver finds a solution +/// outside the element. +std::tuple>> +compute_3D_ray(const dolfinx::mesh::Mesh& mesh, + const xt::xtensor_fixed>& point, + const xt::xtensor_fixed>& tangents, + const std::vector>& cells, + const int max_iter = 25, const double tol = 1e-8); + +} // namespace dolfinx_contact \ No newline at end of file diff --git a/python/dolfinx_contact/wrappers.cpp b/python/dolfinx_contact/wrappers.cpp index 4cc323c6..a970d692 100644 --- a/python/dolfinx_contact/wrappers.cpp +++ b/python/dolfinx_contact/wrappers.cpp @@ -11,6 +11,7 @@ #include #include #include +#include #include #include #include @@ -346,4 +347,41 @@ PYBIND11_MODULE(cpp, m) active_entities); return output; }); + + m.def( + "compute_3D_ray", + [](const dolfinx::mesh::Mesh& mesh, + py::array_t& point, + py::array_t& tangents, + py::array_t& cells, int max_iter, + double tol) + { + // FIXME: How to avoid copy here + auto ents = cells.unchecked(); + const std::size_t shape_0 = cells.shape(0); + std::vector> facets; + facets.reserve(shape_0); + for (std::size_t i = 0; i < shape_0; i++) + facets.emplace_back(ents(i, 0), ents(i, 1)); + + std::array s_p = {(std::size_t)point.shape(0)}; + auto _point + = xt::adapt(point.data(), point.size(), xt::no_ownership(), s_p); + std::array s_t; + std::copy_n(tangents.shape(), 2, s_t.begin()); + auto _tangents = xt::adapt(tangents.data(), tangents.size(), + xt::no_ownership(), s_t); + std::tuple>> + output = dolfinx_contact::compute_3D_ray(mesh, _point, _tangents, + facets, max_iter, tol); + int status = std::get<0>(output); + xt::xtensor x = std::get<2>(output); + std::int32_t idx = std::get<1>(output); + + return py::make_tuple(status, idx, + dolfinx_wrappers::xt_as_pyarray(std::move(x))); + }, + py::arg("mesh"), py::arg("point"), py::arg("tangents"), py::arg("cells"), + py::arg("max_iter") = 25, py::arg("tol") = 1e-8); } diff --git a/python/tests/test_raytracing.py b/python/tests/test_raytracing.py new file mode 100644 index 00000000..e6faab1c --- /dev/null +++ b/python/tests/test_raytracing.py @@ -0,0 +1,77 @@ +# Copyright (C) 2022 Jørgen S. Dokken +# +# SPDX-License-Identifier: MIT +# +# This test check that the ray-tracing routines give the same result as the closest point projection +# when a solution is found (It is not guaranteed that a ray hits the mesh on all processes in parallel) + +import dolfinx_contact +from mpi4py import MPI +import dolfinx +import numpy as np +import pytest + + +@pytest.mark.parametrize("cell_type", [dolfinx.mesh.CellType.hexahedron, dolfinx.mesh.CellType.tetrahedron]) +def test_raytracing(cell_type): + origin = [0.51, 0.33, -1] + tangent = [[1, 0, 0], [0, 1, 0]] + + mesh = dolfinx.mesh.create_unit_cube(MPI.COMM_WORLD, 15, 15, 15, cell_type) + tdim = mesh.topology.dim + facets = dolfinx.mesh.locate_entities_boundary(mesh, tdim - 1, lambda x: np.isclose(x[2], 0)) + + integral_pairs = dolfinx_contact.cpp.compute_active_entities(mesh, facets, dolfinx.fem.IntegralType.exterior_facet) + + status, cell_idx, points = dolfinx_contact.cpp.compute_3D_ray( + mesh, origin, tangent, integral_pairs, 10, 1e-6) + + if status > 0: + # Create structures needed for closest point projections + boundary_cells = dolfinx.mesh.compute_incident_entities(mesh, facets, tdim - 1, tdim) + bb_tree = dolfinx.geometry.BoundingBoxTree(mesh, tdim, boundary_cells) + midpoint_tree = dolfinx.cpp.geometry.create_midpoint_tree(mesh, tdim, boundary_cells) + + # Find closest cell using closest point projection + closest_cell = dolfinx.geometry.compute_closest_entity( + bb_tree, midpoint_tree, mesh, np.reshape(origin, (1, 3)))[0] + assert integral_pairs[cell_idx][0] == closest_cell + + # Compute actual distance between cell and point using GJK + cell_dofs = mesh.geometry.dofmap.links(closest_cell) + cell_geometry = np.empty((len(cell_dofs), 3), dtype=np.float64) + for i, dof in enumerate(cell_dofs): + cell_geometry[i, :] = mesh.geometry.x[dof, :] + distance = dolfinx.geometry.compute_distance_gjk(cell_geometry, origin) + assert np.allclose(points[0], origin + distance) + + +@pytest.mark.parametrize("cell_type", [dolfinx.mesh.CellType.hexahedron, dolfinx.mesh.CellType.tetrahedron]) +def test_raytracing_corner(cell_type): + origin = [1.5, 1.5, 1.5] + tangent = [[1 / np.sqrt(2), 0, -1 / np.sqrt(2)], [-1 / np.sqrt(6), 2 / np.sqrt(6), -1 / np.sqrt(6)]] + + mesh = dolfinx.mesh.create_unit_cube(MPI.COMM_WORLD, 13, 11, 12, cell_type) + tdim = mesh.topology.dim + facets = dolfinx.mesh.locate_entities_boundary(mesh, tdim - 1, lambda x: np.isclose(x[2], 1)) + + integral_pairs = dolfinx_contact.cpp.compute_active_entities(mesh, facets, dolfinx.fem.IntegralType.exterior_facet) + status, cell_idx, points = dolfinx_contact.cpp.compute_3D_ray( + mesh, origin, tangent, integral_pairs, 10, 1e-6) + if status > 0: + # Create structures needed for closest point projections + boundary_cells = dolfinx.mesh.compute_incident_entities(mesh, facets, tdim - 1, tdim) + bb_tree = dolfinx.geometry.BoundingBoxTree(mesh, tdim, boundary_cells) + midpoint_tree = dolfinx.cpp.geometry.create_midpoint_tree(mesh, tdim, boundary_cells) + + # Find closest cell using closest point projection + closest_cell = dolfinx.geometry.compute_closest_entity( + bb_tree, midpoint_tree, mesh, np.reshape(origin, (1, 3)))[0] + + # Compute actual distance between cell and point using GJK + cell_dofs = mesh.geometry.dofmap.links(closest_cell) + cell_geometry = np.empty((len(cell_dofs), 3), dtype=np.float64) + for i, dof in enumerate(cell_dofs): + cell_geometry[i, :] = mesh.geometry.x[dof, :] + distance = dolfinx.geometry.compute_distance_gjk(cell_geometry, origin) + assert np.allclose(points[0], origin + distance)