Skip to content

Commit c947c55

Browse files
authored
Merge branch 'main' into garth/generalise-vector-container
2 parents 463e709 + 855585d commit c947c55

File tree

4 files changed

+32
-3
lines changed

4 files changed

+32
-3
lines changed

python/dolfinx/fem/function.py

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -38,7 +38,9 @@ class Constant(ufl.Constant):
3838
]
3939

4040
def __init__(
41-
self, domain, c: typing.Union[np.ndarray, Sequence, np.floating, np.complexfloating]
41+
self,
42+
domain,
43+
c: typing.Union[float, np.floating, complex, np.complexfloating, Sequence, np.ndarray],
4244
):
4345
"""A constant with respect to a domain.
4446

python/dolfinx/wrappers/fem.cpp

Lines changed: 10 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -103,7 +103,16 @@ void declare_function_space(nb::module_& m, std::string type)
103103
std::shared_ptr<const dolfinx::fem::FiniteElement<T>>>,
104104
std::vector<std::shared_ptr<const dolfinx::fem::DofMap>>>(),
105105
nb::arg("mesh"), nb::arg("elements"), nb::arg("dofmaps"))
106-
.def("collapse", &dolfinx::fem::FunctionSpace<T>::collapse)
106+
.def("collapse",
107+
[](const dolfinx::fem::FunctionSpace<T>& self)
108+
-> std::pair<dolfinx::fem::FunctionSpace<T>,
109+
nanobind::ndarray<std::int32_t, nanobind::numpy>>
110+
{
111+
auto&& [collapsed_fs, dofs] = self.collapse();
112+
return {std::move(collapsed_fs),
113+
dolfinx_wrappers::as_nbarray(std::move(dofs),
114+
{dofs.size()})};
115+
})
107116
.def("component", &dolfinx::fem::FunctionSpace<T>::component)
108117
.def("contains", &dolfinx::fem::FunctionSpace<T>::contains,
109118
nb::arg("V"))

python/dolfinx/wrappers/mesh.cpp

Lines changed: 14 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -696,6 +696,20 @@ void mesh(nb::module_& m)
696696
original_cell_indices.data() + original_cell_indices.size());
697697
},
698698
nb::arg("original_cell_indices"))
699+
.def_prop_ro(
700+
"original_cell_indices",
701+
[](const dolfinx::mesh::Topology& self)
702+
{
703+
const std::vector<std::vector<std::int64_t>>& indices
704+
= self.original_cell_index;
705+
std::vector<nb::ndarray<const std::int64_t, nb::numpy>> idx_nb;
706+
for (auto& oci : indices)
707+
{
708+
idx_nb.push_back(nb::ndarray<const std::int64_t, nb::numpy>(
709+
oci.data(), {oci.size()}));
710+
}
711+
return idx_nb;
712+
})
699713
.def("connectivity",
700714
nb::overload_cast<int, int>(&dolfinx::mesh::Topology::connectivity,
701715
nb::const_),

python/test/unit/fem/test_function_space.py

Lines changed: 5 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -170,9 +170,13 @@ def test_collapse(W, V):
170170
Function(W.sub(1))
171171

172172
Ws = [W.sub(i).collapse() for i in range(W.num_sub_spaces)]
173-
for Wi, _ in Ws:
173+
for Wi, dofs in Ws:
174174
assert np.allclose(Wi.dofmap.index_map.ghosts, W.dofmap.index_map.ghosts)
175175

176+
# Number of collapsed dofs in W numbering must agree with the number of dofs
177+
# of the collapsed space
178+
assert Wi.dofmap.index_map.size_local + Wi.dofmap.index_map.num_ghosts == dofs.size
179+
176180
msh = W.mesh
177181
cell_imap = msh.topology.index_map(msh.topology.dim)
178182
num_cells = cell_imap.size_local + cell_imap.num_ghosts

0 commit comments

Comments
 (0)