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
BoxMesh.h
1 // Copyright (C) 2005-2017 Anders Logg and 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 <array>
10 #include <cstddef>
11 #include <dolfinx/graph/AdjacencyList.h>
12 #include <dolfinx/mesh/Mesh.h>
13 #include <dolfinx/mesh/cell_types.h>
14 #include <mpi.h>
15 
16 namespace dolfinx::fem
17 {
18 class CoordinateElement;
19 }
20 
23 {
24 
41 create(MPI_Comm comm, const std::array<std::array<double, 3>, 2>& p,
42  std::array<std::size_t, 3> n, mesh::CellType celltype,
43  const mesh::GhostMode ghost_mode,
44  const mesh::CellPartitionFunction& partitioner
45  = static_cast<graph::AdjacencyList<std::int32_t> (*)(
46  MPI_Comm, int, int, const graph::AdjacencyList<std::int64_t>&,
48 } // namespace dolfinx::generation::BoxMesh
This class provides a static adjacency list data structure. It is commonly used to store directed gra...
Definition: AdjacencyList.h:47
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:53
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
Right cuboid mesh creation.
Definition: BoxMesh.h:23
mesh::Mesh create(MPI_Comm comm, const std::array< std::array< double, 3 >, 2 > &p, std::array< std::size_t, 3 > n, mesh::CellType celltype, const mesh::GhostMode ghost_mode, const mesh::CellPartitionFunction &partitioner=static_cast< graph::AdjacencyList< std::int32_t >(*)(MPI_Comm, int, int, const graph::AdjacencyList< std::int64_t > &, mesh::GhostMode)>(&mesh::partition_cells_graph))
Create a uniform mesh::Mesh over the rectangular prism spanned by the two points p....
Definition: BoxMesh.cpp:203
graph::AdjacencyList< std::int32_t > partition_cells_graph(MPI_Comm comm, int n, int tdim, const graph::AdjacencyList< std::int64_t > &cells, mesh::GhostMode ghost_mode)
Compute destination rank for mesh cells in this rank by applying the default graph partitioner to the...
Definition: utils.cpp:504
std::function< const dolfinx::graph::AdjacencyList< std::int32_t >(MPI_Comm comm, int nparts, int tdim, const dolfinx::graph::AdjacencyList< std::int64_t > &cells, dolfinx::mesh::GhostMode ghost_mode)> CellPartitionFunction
Definition: Mesh.h:40
CellType
Cell type identifier.
Definition: cell_types.h:22
GhostMode
Enum for different partitioning ghost modes.
Definition: Mesh.h:44