Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.3.0/v0.9.0/cpp
DOLFINx  0.3.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::mesh
15 {
16 class Mesh;
17 }
18 
19 namespace dolfinx::geometry
20 {
21 class BoundingBoxTree;
22 
29 BoundingBoxTree
30 create_midpoint_tree(const mesh::Mesh& mesh, int tdim,
31  const xtl::span<const std::int32_t>& entity_indices);
32 
38 std::vector<std::array<int, 2>>
39 compute_collisions(const BoundingBoxTree& tree0, const BoundingBoxTree& tree1);
40 
45 std::vector<int> compute_collisions(const BoundingBoxTree& tree,
46  const std::array<double, 3>& p);
47 
56 std::pair<std::int32_t, double>
57 compute_closest_entity(const BoundingBoxTree& tree,
58  const std::array<double, 3>& p, const mesh::Mesh& mesh,
59  double R = -1);
60 
63 double
64 compute_squared_distance_bbox(const std::array<std::array<double, 3>, 2>& b,
65  const std::array<double, 3>& x);
66 
80 double squared_distance(const mesh::Mesh& mesh, int dim, std::int32_t index,
81  const std::array<double, 3>& p);
82 
91 std::vector<std::int32_t>
93  const xtl::span<const std::int32_t>& candidate_cells,
94  const std::array<double, 3>& p, int n);
95 } // namespace dolfinx::geometry
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:53
Geometry data structures and algorithms.
Definition: BoundingBoxTree.h:20
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
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
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
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
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
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
Mesh data structures and algorithms on meshes.
Definition: DirichletBC.h:20