DOLFIN-X
DOLFIN-X C++ interface
|
10 #include <dolfinx/graph/AdjacencyList.h>
16 class ElementDofLayout;
33 graph::AdjacencyList<std::int64_t>
35 const graph::AdjacencyList<std::int64_t>& cells);
39 const Eigen::Ref<const Eigen::ArrayXi>& entities,
44 const Eigen::Ref<const Eigen::ArrayXi>& entities,
48 Eigen::ArrayXd
h(
const Mesh& mesh,
49 const Eigen::Ref<const Eigen::ArrayXi>& entities,
int dim);
52 Eigen::ArrayXd
inradius(
const Mesh& mesh,
53 const Eigen::Ref<const Eigen::ArrayXi>& entities);
57 const Eigen::Ref<const Eigen::ArrayXi>& entities);
60 Eigen::Array<double, Eigen::Dynamic, 3, Eigen::RowMajor>
64 Eigen::Array<double, Eigen::Dynamic, 3, Eigen::RowMajor>
midpoints(
65 const mesh::Mesh& mesh,
int dim,
66 const Eigen::Ref<
const Eigen::Array<int, Eigen::Dynamic, 1>>& entities);
79 const mesh::Mesh& mesh,
const int dim,
80 const std::function<Eigen::Array<bool, Eigen::Dynamic, 1>(
81 const Eigen::Ref<
const Eigen::Array<
double, 3, Eigen::Dynamic,
82 Eigen::RowMajor>>&)>& marker);
106 const mesh::Mesh& mesh,
const int dim,
107 const std::function<Eigen::Array<bool, Eigen::Dynamic, 1>(
108 const Eigen::Ref<
const Eigen::Array<
double, 3, Eigen::Dynamic,
109 Eigen::RowMajor>>&)>& marker);
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:321
CellType
Cell type identifier.
Definition: cell_types.h:23
Eigen::ArrayXd circumradius(const Mesh &mesh, const Eigen::Ref< const Eigen::ArrayXi > &entities, int dim)
Compute circumradius of mesh entities.
Definition: utils.cpp:393
Eigen::Array< double, Eigen::Dynamic, 3, Eigen::RowMajor > cell_normals(const Mesh &mesh, int dim)
Compute normal to given cell (viewed as embedded in 3D)
Definition: utils.cpp:462
Eigen::ArrayXd volume_entities(const Mesh &mesh, const Eigen::Ref< const Eigen::ArrayXi > &entities, int dim)
Compute (generalized) volume of mesh entities of given dimension.
Definition: utils.cpp:350
Eigen::ArrayXd h(const Mesh &mesh, const Eigen::Ref< const Eigen::ArrayXi > &entities, int dim)
Compute greatest distance between any two vertices.
Definition: utils.cpp:356
Eigen::Array< double, Eigen::Dynamic, 3, Eigen::RowMajor > midpoints(const mesh::Mesh &mesh, int dim, const Eigen::Ref< const Eigen::Array< int, Eigen::Dynamic, 1 >> &entities)
Compute midpoints or mesh entities of a given dimension.
Definition: utils.cpp:538
Eigen::Array< std::int32_t, Eigen::Dynamic, 1 > locate_entities_boundary(const mesh::Mesh &mesh, const int dim, const std::function< Eigen::Array< bool, Eigen::Dynamic, 1 >(const Eigen::Ref< const Eigen::Array< double, 3, Eigen::Dynamic, Eigen::RowMajor >> &)> &marker)
Compute indicies of all mesh entities that are attached to an owned boundary facet and evaluate to tr...
Definition: utils.cpp:670
Eigen::ArrayXd radius_ratio(const Mesh &mesh, const Eigen::Ref< const Eigen::ArrayXi > &entities)
Compute dim*inradius/circumradius for given cells.
Definition: utils.cpp:451
Eigen::ArrayXd inradius(const Mesh &mesh, const Eigen::Ref< const Eigen::ArrayXi > &entities)
Compute inradius of cells.
Definition: utils.cpp:399
Eigen::Array< std::int32_t, Eigen::Dynamic, 1 > locate_entities(const mesh::Mesh &mesh, const int dim, const std::function< Eigen::Array< bool, Eigen::Dynamic, 1 >(const Eigen::Ref< const Eigen::Array< double, 3, Eigen::Dynamic, Eigen::RowMajor >> &)> &marker)
Compute indicies of all mesh entities that evaluate to true for the provided geometric marking functi...
Definition: utils.cpp:601