DOLFIN-X
DOLFIN-X C++ interface
|
10 #include <dolfinx/graph/AdjacencyList.h>
16 class ElementDofLayout;
34 graph::AdjacencyList<std::int64_t>
36 const graph::AdjacencyList<std::int64_t>& cells);
40 const Eigen::Ref<const Eigen::ArrayXi>& entities,
45 const Eigen::Ref<const Eigen::ArrayXi>& entities,
49 Eigen::ArrayXd
h(
const Mesh& mesh,
50 const Eigen::Ref<const Eigen::ArrayXi>& entities,
int dim);
53 Eigen::ArrayXd
inradius(
const Mesh& mesh,
54 const Eigen::Ref<const Eigen::ArrayXi>& entities);
58 const Eigen::Ref<const Eigen::ArrayXi>& entities);
61 Eigen::Array<double, Eigen::Dynamic, 3, Eigen::RowMajor>
65 Eigen::Vector3d
normal(
const MeshEntity& cell,
int facet_local);
68 Eigen::Array<double, Eigen::Dynamic, 3, Eigen::RowMajor>
midpoints(
69 const mesh::Mesh& mesh,
int dim,
70 const Eigen::Ref<
const Eigen::Array<int, Eigen::Dynamic, 1>>& entities);
85 const mesh::Mesh& mesh,
const int dim,
86 const std::function<Eigen::Array<bool, Eigen::Dynamic, 1>(
87 const Eigen::Ref<
const Eigen::Array<
double, 3, Eigen::Dynamic,
88 Eigen::RowMajor>>&)>& marker,
89 const bool boundary_only);
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:22
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< std::int32_t, Eigen::Dynamic, 1 > locate_entities_geometrical(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, const bool boundary_only)
Compute indicies (local to the process) of all mesh entities that evaluate to true for the provided m...
Definition: utils.cpp:728
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:665
Eigen::Vector3d normal(const MeshEntity &cell, int facet_local)
Compute of given facet with respect to the cell.
Definition: utils.cpp:541
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