DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
utils.h File Reference

Functions supporting mesh operations. More...

#include "Mesh.h"
#include "Topology.h"
#include "graphbuild.h"
#include <algorithm>
#include <basix/mdspan.hpp>
#include <concepts>
#include <dolfinx/graph/AdjacencyList.h>
#include <dolfinx/graph/ordering.h>
#include <dolfinx/graph/partition.h>
#include <functional>
#include <mpi.h>
#include <numeric>
#include <span>
Include dependency graph for utils.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Namespaces

namespace  dolfinx
 Top-level namespace.
 
namespace  dolfinx::fem
 Finite element method functionality.
 
namespace  dolfinx::mesh
 Mesh data structures and algorithms on meshes.
 

Concepts

concept  dolfinx::mesh::MarkerFn
 Requirements on function for geometry marking.
 

Typedefs

using CellPartitionFunction
 Signature for the cell partitioning function. Function that implement this interface compute the destination rank for cells currently on this rank.
 

Enumerations

enum class  GhostMode : int { none , shared_facet , shared_vertex }
 Enum for different partitioning ghost modes.
 

Functions

template<typename T >
void reorder_list (std::span< T > list, std::span< const std::int32_t > nodemap)
 Re-order the nodes of a fixed-degree adjacency list.
 
template<std::floating_point T>
std::tuple< std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > compute_vertex_coords_boundary (const mesh::Mesh< T > &mesh, int dim, std::span< const std::int32_t > facets)
 Compute the coordinates of 'vertices' for entities of a given dimension that are attached to specified facets.
 
std::vector< std::int32_t > exterior_facet_indices (const Topology &topology)
 Compute the indices of all exterior facets that are owned by the caller.
 
std::vector< std::int64_t > extract_topology (CellType cell_type, const fem::ElementDofLayout &layout, std::span< const std::int64_t > cells)
 Extract topology from cell data, i.e. extract cell vertices.
 
template<std::floating_point T>
std::vector< T > h (const Mesh< T > &mesh, std::span< const std::int32_t > entities, int dim)
 Compute greatest distance between any two vertices of the mesh entities (h).
 
template<std::floating_point T>
std::vector< T > cell_normals (const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities)
 Compute normal to given cell (viewed as embedded in 3D).
 
template<std::floating_point T>
std::vector< T > compute_midpoints (const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities)
 Compute the midpoints for mesh entities of a given dimension.
 
template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 2 > > compute_vertex_coords (const mesh::Mesh< T > &mesh)
 The coordinates for all 'vertices' in the mesh.
 
template<std::floating_point T, MarkerFn< T > U>
std::vector< std::int32_t > locate_entities (const Mesh< T > &mesh, int dim, U marker)
 Compute indices of all mesh entities that evaluate to true for the provided geometric marking function.
 
template<std::floating_point T, MarkerFn< T > U>
std::vector< std::int32_t > locate_entities_boundary (const Mesh< T > &mesh, int dim, U marker)
 Compute indices of all mesh entities that are attached to an owned boundary facet and evaluate to true for the provided geometric marking function.
 
template<std::floating_point T>
std::vector< std::int32_t > entities_to_geometry (const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities, bool permute=false)
 Compute the geometry degrees of freedom associated with the closure of a given set of cell entities.
 
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 on this rank by applying the default graph partitioner to the dual graph of the mesh.
 
std::vector< std::int32_t > compute_incident_entities (const Topology &topology, std::span< const std::int32_t > entities, int d0, int d1)
 Compute incident entities.
 
template<typename U >
Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh (MPI_Comm comm, MPI_Comm commt, std::vector< std::span< const std::int64_t > > cells, const std::vector< fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > > &elements, MPI_Comm commg, const U &x, std::array< std::size_t, 2 > xshape, const CellPartitionFunction &partitioner)
 Create a distributed mesh::Mesh from mesh data and using the provided graph partitioning function for determining the parallel distribution of the mesh.
 
template<typename U >
Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh (MPI_Comm comm, MPI_Comm commt, std::span< const std::int64_t > cells, const fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > &element, MPI_Comm commg, const U &x, std::array< std::size_t, 2 > xshape, const CellPartitionFunction &partitioner)
 Create a distributed mesh with a single cell type from mesh data and using a provided graph partitioning function for determining the parallel distribution of the mesh.
 
template<typename U >
Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh (MPI_Comm comm, std::span< const std::int64_t > cells, const fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > &elements, const U &x, std::array< std::size_t, 2 > xshape, GhostMode ghost_mode)
 Create a distributed mesh from mesh data using the default graph partitioner to determine the parallel distribution of the mesh.
 
template<std::floating_point T>
std::pair< Geometry< T >, std::vector< int32_t > > create_subgeometry (const Mesh< T > &mesh, int dim, std::span< const std::int32_t > subentity_to_entity)
 Create a sub-geometry from a mesh and a subset of mesh entities to be included.
 
template<std::floating_point T>
std::tuple< Mesh< T >, std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< std::int32_t > > create_submesh (const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities)
 Create a new mesh consisting of a subset of entities in a mesh.
 

Detailed Description

Functions supporting mesh operations.

Function Documentation

◆ compute_vertex_coords()

template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 2 > > compute_vertex_coords ( const mesh::Mesh< T > & mesh)

The coordinates for all 'vertices' in the mesh.

Parameters
[in]meshMesh to compute the vertex coordinates for.
Returns
The vertex coordinates. The shape is (3, num_vertices) and the jth column hold the coordinates of vertex j.

◆ compute_vertex_coords_boundary()

template<std::floating_point T>
std::tuple< std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > compute_vertex_coords_boundary ( const mesh::Mesh< T > & mesh,
int dim,
std::span< const std::int32_t > facets )

Compute the coordinates of 'vertices' for entities of a given dimension that are attached to specified facets.

Precondition
The provided facets must be on the boundary of the mesh.
Parameters
[in]meshMesh to compute the vertex coordinates for.
[in]dimTopological dimension of the entities.
[in]facetsList of facets (must be on the mesh boundary).
Returns
(0) Entities attached to the boundary facets (sorted), (1) vertex coordinates (shape is (3, num_vertices)) and (2) map from vertex in the full mesh to the position in the vertex coordinates array (set to -1 if vertex in full mesh is not in the coordinate array).

◆ reorder_list()

template<typename T >
void reorder_list ( std::span< T > list,
std::span< const std::int32_t > nodemap )

Re-order the nodes of a fixed-degree adjacency list.

Parameters
[in,out]listFixed-degree adjacency list stored row-major. Degree is equal to list.size() / nodemap.size().
[in]nodemapMap from old to new index, i.e. for an old index i the new index is nodemap[i].