Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d0/d47/mesh_2utils_8h_source.html
DOLFINx 0.6.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
utils.h
1// Copyright (C) 2019-2020 Garth N. Wells
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 <basix/mdspan.hpp>
10#include <dolfinx/graph/AdjacencyList.h>
11#include <dolfinx/graph/partition.h>
12#include <functional>
13#include <mpi.h>
14#include <span>
15
16namespace dolfinx::fem
17{
18class ElementDofLayout;
19}
20
21namespace dolfinx::mesh
22{
23enum class CellType;
24class Mesh;
25class Topology;
26
28enum class GhostMode : int
29{
30 none,
31 shared_facet,
32 shared_vertex
33};
34
35// See https://github.com/doxygen/doxygen/issues/9552
38/*
52*/
53using CellPartitionFunction = std::function<graph::AdjacencyList<std::int32_t>(
54 MPI_Comm comm, int nparts, int tdim,
56
67extract_topology(const CellType& cell_type, const fem::ElementDofLayout& layout,
69
78std::vector<double> h(const Mesh& mesh, std::span<const std::int32_t> entities,
79 int dim);
80
84std::vector<double> cell_normals(const Mesh& mesh, int dim,
85 std::span<const std::int32_t> entities);
86
90std::vector<double> compute_midpoints(const Mesh& mesh, int dim,
91 std::span<const std::int32_t> entities);
92
103std::vector<std::int32_t> locate_entities(
104 const Mesh& mesh, int dim,
105 const std::function<std::vector<std::int8_t>(
106 std::experimental::mdspan<
107 const double,
108 std::experimental::extents<
109 std::size_t, 3, std::experimental::dynamic_extent>>)>& marker);
110
131std::vector<std::int32_t> locate_entities_boundary(
132 const Mesh& mesh, int dim,
133 const std::function<std::vector<std::int8_t>(
134 std::experimental::mdspan<
135 const double,
136 std::experimental::extents<
137 std::size_t, 3, std::experimental::dynamic_extent>>)>& marker);
138
155std::vector<std::int32_t>
156entities_to_geometry(const Mesh& mesh, int dim,
157 std::span<const std::int32_t> entities, bool orient);
158
170std::vector<std::int32_t> exterior_facet_indices(const Topology& topology);
171
177 = mesh::GhostMode::none,
178 const graph::partition_fn& partfn
180
188std::vector<std::int32_t> compute_incident_entities(
189 const Mesh& mesh, std::span<const std::int32_t> entities, int d0, int d1);
190
191} // namespace dolfinx::mesh
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh ent...
Definition: ElementDofLayout.h:31
This class provides a static adjacency list data structure. It is commonly used to store directed gra...
Definition: AdjacencyList.h:27
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:34
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relat...
Definition: Topology.h:43
Finite element method functionality.
Definition: assemble_matrix_impl.h:25
std::function< graph::AdjacencyList< std::int32_t >(MPI_Comm, int, const AdjacencyList< std::int64_t > &, bool)> partition_fn
Signature of functions for computing the parallel partitioning of a distributed graph.
Definition: partition.h:36
AdjacencyList< std::int32_t > partition_graph(MPI_Comm comm, int nparts, const AdjacencyList< std::int64_t > &local_graph, bool ghosting)
Partition graph across processes using the default graph partitioner.
Definition: partition.cpp:21
Mesh data structures and algorithms on meshes.
Definition: DofMap.h:31
GhostMode
Enum for different partitioning ghost modes.
Definition: utils.h:29
std::vector< std::int32_t > compute_incident_entities(const Mesh &mesh, std::span< const std::int32_t > entities, int d0, int d1)
Compute incident indices.
Definition: utils.cpp:630
std::function< graph::AdjacencyList< std::int32_t >(MPI_Comm comm, int nparts, int tdim, const graph::AdjacencyList< std::int64_t > &cells)> CellPartitionFunction
Signature for the cell partitioning function. The function should compute the destination rank for ce...
Definition: utils.h:55
std::vector< double > compute_midpoints(const Mesh &mesh, int dim, std::span< const std::int32_t > entities)
Compute the midpoints for mesh entities of a given dimension.
Definition: utils.cpp:357
std::vector< std::int32_t > locate_entities_boundary(const Mesh &mesh, int dim, const std::function< std::vector< std::int8_t >(std::experimental::mdspan< const double, std::experimental::extents< std::size_t, 3, std::experimental::dynamic_extent > >)> &marker)
Compute indices of all mesh entities that are attached to an owned boundary facet and evaluate to tru...
Definition: utils.cpp:436
std::vector< std::int32_t > locate_entities(const Mesh &mesh, int dim, const std::function< std::vector< std::int8_t >(std::experimental::mdspan< const double, std::experimental::extents< std::size_t, 3, std::experimental::dynamic_extent > >)> &marker)
Compute indices of all mesh entities that evaluate to true for the provided geometric marking functio...
Definition: utils.cpp:388
std::vector< std::int32_t > entities_to_geometry(const Mesh &mesh, int dim, std::span< const std::int32_t > entities, bool orient)
Determine the indices in the geometry data for each vertex of the given mesh entities.
Definition: utils.cpp:494
graph::AdjacencyList< std::int64_t > extract_topology(const CellType &cell_type, const fem::ElementDofLayout &layout, const graph::AdjacencyList< std::int64_t > &cells)
Extract topology from cell data, i.e. extract cell vertices.
Definition: utils.cpp:163
std::vector< std::int32_t > exterior_facet_indices(const Topology &topology)
Compute the indices of all exterior facets that are owned by the caller.
Definition: utils.cpp:580
CellType
Cell type identifier.
Definition: cell_types.h:22
std::vector< double > cell_normals(const Mesh &mesh, int dim, std::span< const std::int32_t > entities)
Compute normal to given cell (viewed as embedded in 3D)
Definition: utils.cpp:240
std::vector< double > h(const Mesh &mesh, std::span< const std::int32_t > entities, int dim)
Compute greatest distance between any two vertices of the mesh entities (h).
Definition: utils.cpp:190
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph)
Create a function that computes destination rank for mesh cells in this rank by applying the default ...
Definition: utils.cpp:609