Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.0
DOLFINx C++ interface
utils.h
1 // Copyright (C) 2019-2021 Garth N. Wells and Jørgen S. Dokken
2 //
3 // This file is part of DOLFINx (https://www.fenicsproject.org)
4 //
5 // SPDX-License-Identifier: LGPL-3.0-or-later
6 
7 #pragma once
8 
9 #include <array>
10 #include <utility>
11 #include <vector>
12 #include <xtl/xspan.hpp>
13 
14 namespace dolfinx
15 {
16 namespace mesh
17 {
18 class Mesh;
19 } // namespace mesh
20 
21 namespace geometry
22 {
23 class BoundingBoxTree;
24 
31 BoundingBoxTree
32 create_midpoint_tree(const mesh::Mesh& mesh, int tdim,
33  const xtl::span<const std::int32_t>& entity_indices);
34 
40 std::vector<std::array<int, 2>>
41 compute_collisions(const BoundingBoxTree& tree0, const BoundingBoxTree& tree1);
42 
47 std::vector<int> compute_collisions(const BoundingBoxTree& tree,
48  const std::array<double, 3>& p);
49 
58 std::pair<std::int32_t, double>
59 compute_closest_entity(const BoundingBoxTree& tree,
60  const std::array<double, 3>& p, const mesh::Mesh& mesh,
61  double R = -1);
62 
65 double
66 compute_squared_distance_bbox(const std::array<std::array<double, 3>, 2>& b,
67  const std::array<double, 3>& x);
68 
82 double squared_distance(const mesh::Mesh& mesh, int dim, std::int32_t index,
83  const std::array<double, 3>& p);
84 
93 std::vector<std::int32_t>
95  const xtl::span<const std::int32_t>& candidate_cells,
96  const std::array<double, 3>& p, int n);
97 } // namespace geometry
98 } // namespace dolfinx
dolfinx::geometry::squared_distance
double squared_distance(const mesh::Mesh &mesh, int dim, std::int32_t index, const std::array< double, 3 > &p)
Compute squared distance from a given point to the nearest point on a cell (only first order convex c...
Definition: utils.cpp:308
dolfinx::geometry::select_colliding_cells
std::vector< std::int32_t > select_colliding_cells(const dolfinx::mesh::Mesh &mesh, const xtl::span< const std::int32_t > &candidate_cells, const std::array< double, 3 > &p, int n)
From the given Mesh, select up to n cells (local to process) from the list which actually collide wit...
Definition: utils.cpp:364
dolfinx::mesh::Mesh
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:55
dolfinx::geometry::compute_closest_entity
std::pair< std::int32_t, double > compute_closest_entity(const BoundingBoxTree &tree, const std::array< double, 3 > &p, const mesh::Mesh &mesh, double R=-1)
Compute closest mesh entity (local to process) for the topological distance of the bounding box tree ...
Definition: utils.cpp:271
dolfinx::geometry::create_midpoint_tree
BoundingBoxTree create_midpoint_tree(const mesh::Mesh &mesh, int tdim, const xtl::span< const std::int32_t > &entity_indices)
Create a bounding box tree for a subset of entities (local to process) based on the entity midpoints.
Definition: utils.cpp:209
dolfinx::geometry::compute_collisions
std::vector< std::array< int, 2 > > compute_collisions(const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
Compute all collisions between two BoundingBoxTrees (local to process).
Definition: utils.cpp:230
dolfinx::geometry::compute_squared_distance_bbox
double compute_squared_distance_bbox(const std::array< std::array< double, 3 >, 2 > &b, const std::array< double, 3 > &x)
Compute squared distance between point and bounding box wih index "node". Returns zero if point is in...
Definition: utils.cpp:254