Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Nanobind + Pyproject.toml #91

Merged
merged 30 commits into from
Nov 20, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
Show all changes
30 commits
Select commit Hold shift + click to select a range
c5953ac
Start reworking mpc with nanobind. Some tests pass. Currently debuggi…
jorgensd Nov 18, 2023
8438260
Fixes of implicit conversion of integers to garbage
jorgensd Nov 18, 2023
58d7a03
Various conversion for nanobind
jorgensd Nov 18, 2023
6fe20bf
Remove consts
jorgensd Nov 18, 2023
c2debf5
Typing
jorgensd Nov 18, 2023
9726797
Deactivate nonlinear test for now
jorgensd Nov 18, 2023
348aee0
Update CI
jorgensd Nov 18, 2023
8a120d4
Turn of ufl import
jorgensd Nov 18, 2023
4f210de
Fix lifting default arg
jorgensd Nov 19, 2023
ce45a24
flake8
jorgensd Nov 19, 2023
5e893bf
Remove muddling build types
jorgensd Nov 19, 2023
5132a7d
Flake8 + mypy
jorgensd Nov 19, 2023
f20a914
Mypy
jorgensd Nov 19, 2023
67c05e3
Lifting type fixes
jorgensd Nov 19, 2023
ecbf72a
Switch if order
jorgensd Nov 19, 2023
271d598
mypy
jorgensd Nov 19, 2023
f23c729
Renable tests due to fixes in https://github.com/FEniCS/dolfinx/pull/…
jorgensd Nov 19, 2023
baf3f79
Fix assemble_vector
jorgensd Nov 19, 2023
8ce3d79
update action
jorgensd Nov 19, 2023
984b34e
More workflows
jorgensd Nov 19, 2023
138d0f6
FIx connectivity
jorgensd Nov 20, 2023
005b8eb
Fix axtion and tolerance in periodic constraint
jorgensd Nov 20, 2023
bed25e8
Tweak action
jorgensd Nov 20, 2023
76517a8
More tweaks
jorgensd Nov 20, 2023
e0f91b7
Various fixes
jorgensd Nov 20, 2023
61b690f
More demo fixes
jorgensd Nov 20, 2023
3c48670
Mypy
jorgensd Nov 20, 2023
c405ff0
Fix edge example
jorgensd Nov 20, 2023
4b9c660
Remove reshape
jorgensd Nov 20, 2023
facbe47
Arange for cells in midpoints
jorgensd Nov 20, 2023
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
4 changes: 3 additions & 1 deletion .github/actions/install-dolfinx/action.yaml
Original file line number Diff line number Diff line change
Expand Up @@ -70,4 +70,6 @@ runs:

- name: Build Python interface (dolfinx)
shell: bash -el {0}
run: PETSC_DIR=${{ inputs.petsc_dir }} PETSC_ARCH=${{ inputs.petsc_arch }} BUILD_TYPE="Release" python3 -m pip -v install ./dolfinx/python/
run: |
python3 -m pip install -r ./dolfinx/python/build-requirements.txt
PETSC_DIR=${{ inputs.petsc_dir }} PETSC_ARCH=${{ inputs.petsc_arch }} python3 -m pip -v install --config-settings=cmake.build-type="Release" --no-build-isolation ./dolfinx/python
2 changes: 1 addition & 1 deletion .github/workflows/build_docs.yml
Original file line number Diff line number Diff line change
Expand Up @@ -36,7 +36,7 @@ jobs:
cmake --install build-dir

- name: Install DOLFINx-MPC (Python)
run: python3 -m pip -v install python/[docs]
run: python3 -m pip -v install --config-settings=cmake.build-type="Release" --no-build-isolation ./python/[docs]

- name: Build docs
run: jupyter book build .
Expand Down
2 changes: 1 addition & 1 deletion .github/workflows/sonarcloud.yml
Original file line number Diff line number Diff line change
Expand Up @@ -81,7 +81,7 @@ jobs:
cmake --install build-dir

- name: Install DOLFINx-MPC (Python)
run: python3 -m pip -v install python/
run: python3 -m pip -v install --config-settings=cmake.build-type="Release" --no-build-isolation ./python

- name: Run sonar-scanner
env:
Expand Down
2 changes: 1 addition & 1 deletion .github/workflows/test_mpc.yml
Original file line number Diff line number Diff line change
Expand Up @@ -99,7 +99,7 @@ jobs:
cmake --install build-dir

- name: Install DOLFINx-MPC (Python)
run: CXX_FLAGS="${MPC_CMAKE_CXX_FLAGS}" python3 -m pip -v install -e python/[test]
run: python3 -m pip -v install --config-settings=cmake.build-type=${MPC_BUILD_MODE} --no-build-isolation -e python/[test]


- name: Run tests (serial)
Expand Down
7 changes: 6 additions & 1 deletion Changelog.md
Original file line number Diff line number Diff line change
@@ -1,6 +1,11 @@
# Changelog

## Main
- **API**
- Various shared pointers in C++ interface is changed to const references
- Multipoint-constraint now accept `std::span` instead of vectors
- Now using [nanobind](https://github.com/wjakob/nanobind) for Python bindings
- Switch to `pyproject.toml`, **see installation notes** for updated instructions
- **DOLFINx API-changes**
- `dolfinx.fem.FunctionSpaceBase` replaced by `dolfinx.fem.FunctionSpace`
- `ufl.FiniteElement` and `ufl.VectorElement` is replaced by `basix.ufl.element`
Expand All @@ -10,7 +15,7 @@
- **API**:
- Change input of `dolfinx_mpc.MultiPointConstraint.homogenize` and `dolfinx_mpc.backsubstitution` to `dolfinx.fem.Function` instead of `PETSc.Vec`.
- **New feature**: Add support for more floating types (float32, float64, complex64, complex128). The floating type of a MPC is related to the mesh geometry.
- This resulted in a minor refactoring of the pybindings, meaning that tte class `dolfinx_mpc.cpp.mpc.MultiPointConstraint` is replaced by `dolfinx_mpc.cpp.mpc.MultiPointConstraint_{dtype}`
- This resulted in a minor refactoring of the pybindings, meaning that the class `dolfinx_mpc.cpp.mpc.MultiPointConstraint` is replaced by `dolfinx_mpc.cpp.mpc.MultiPointConstraint_{dtype}`
- Casting scalar-type with `dolfinx.default_scalar_type` instead of `PETSc.ScalarType`
- Remove usage of `VectorFunctionSpace`. Use blocked basix element instead.
- **DOLFINX API-changes**:
Expand Down
2 changes: 1 addition & 1 deletion README.md
Original file line number Diff line number Diff line change
Expand Up @@ -61,5 +61,5 @@ To install the `dolfinx_mpc`-library run the following code from this directory:
```bash
cmake -G Ninja -DCMAKE_BUILD_TYPE=Release -B build-dir cpp/
ninja -j3 install -C build-dir
python3 -m pip install python/. --upgrade
python3 -m pip -v install --config-settings=cmake.build-type="Release" --no-build-isolation ./python -U
```
88 changes: 43 additions & 45 deletions cpp/ContactConstraint.h
Original file line number Diff line number Diff line change
Expand Up @@ -63,17 +63,16 @@ dolfinx_mpc::mpc_data<T> compute_master_contributions(
U, MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
std::size_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent, 3>>
normals,
std::shared_ptr<const dolfinx::fem::FunctionSpace<U>> V,
const dolfinx::fem::FunctionSpace<U>& V,
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
tabulated_basis_values)
{
const double tol = 1e-6;
auto mesh = V->mesh();
const std::int32_t block_size = V->dofmap()->index_map_bs();
std::shared_ptr<const dolfinx::common::IndexMap> imap
= V->dofmap()->index_map;
const int bs = V->dofmap()->index_map_bs();
auto mesh = V.mesh();
const std::int32_t block_size = V.dofmap()->index_map_bs();
std::shared_ptr<const dolfinx::common::IndexMap> imap = V.dofmap()->index_map;
const int bs = V.dofmap()->index_map_bs();
const std::int32_t size_local = imap->size_local();
MPI_Comm comm = mesh->comm();
int rank = -1;
Expand All @@ -89,7 +88,7 @@ dolfinx_mpc::mpc_data<T> compute_master_contributions(
{
if (const std::int32_t cell = local_colliding_cell[i]; cell != -1)
{
auto cell_blocks = V->dofmap()->cell_dofs(cell);
auto cell_blocks = V.dofmap()->cell_dofs(cell);
for (std::size_t j = 0; j < cell_blocks.size(); ++j)
{
for (int b = 0; b < bs; b++)
Expand Down Expand Up @@ -124,7 +123,7 @@ dolfinx_mpc::mpc_data<T> compute_master_contributions(
{
if (const std::int32_t cell = local_colliding_cell[i]; cell != -1)
{
auto cell_blocks = V->dofmap()->cell_dofs(cell);
auto cell_blocks = V.dofmap()->cell_dofs(cell);
global_blocks.resize(cell_blocks.size());
imap->local_to_global(cell_blocks, global_blocks);
// Compute coefficients for each master
Expand Down Expand Up @@ -169,7 +168,7 @@ dolfinx_mpc::mpc_data<T> compute_master_contributions(
/// tagged with the marker
template <std::floating_point U>
std::vector<std::int32_t>
locate_slave_dofs(std::shared_ptr<const dolfinx::fem::FunctionSpace<U>> V,
locate_slave_dofs(const dolfinx::fem::FunctionSpace<U>& V,
const dolfinx::mesh::MeshTags<std::int32_t>& meshtags,
std::int32_t slave_marker)
{
Expand All @@ -184,13 +183,12 @@ locate_slave_dofs(std::shared_ptr<const dolfinx::fem::FunctionSpace<U>> V,
// Find all dofs on slave facets
// NOTE: Assumption that we are only working with vector spaces, which is
// ordered as xyz,xyzgeometry
auto [V0, map] = V->sub({0})->collapse();
auto [V0, map] = V.sub({0})->collapse();

std::array<std::vector<std::int32_t>, 2> slave_dofs
= dolfinx::fem::locate_dofs_topological(
*V->mesh()->topology_mutable(),
{*V->sub({0})->dofmap(), *V0.dofmap()}, edim,
std::span(slave_facets));
*V.mesh()->topology_mutable(), {*V.sub({0})->dofmap(), *V0.dofmap()},
edim, std::span(slave_facets));
return slave_dofs[0];
}

Expand Down Expand Up @@ -347,27 +345,27 @@ namespace dolfinx_mpc
/// collision
template <typename T, std::floating_point U>
mpc_data<T> create_contact_slip_condition(
std::shared_ptr<dolfinx::fem::FunctionSpace<U>> V,
const dolfinx::fem::FunctionSpace<U>& V,
const dolfinx::mesh::MeshTags<std::int32_t>& meshtags,
std::int32_t slave_marker, std::int32_t master_marker,
std::shared_ptr<dolfinx::fem::Function<T, U>> nh, const U eps2 = 1e-20)
const dolfinx::fem::Function<T, U>& nh, const U eps2 = 1e-20)
{

dolfinx::common::Timer timer("~MPC: Create slip constraint");
std::shared_ptr<const mesh::Mesh<U>> mesh = V->mesh();
std::shared_ptr<const mesh::Mesh<U>> mesh = V.mesh();
MPI_Comm comm = mesh->comm();
int rank = -1;
MPI_Comm_rank(comm, &rank);

// Extract some const information from function-space
const std::shared_ptr<const dolfinx::common::IndexMap> imap
= V->dofmap()->index_map;
= V.dofmap()->index_map;
assert(mesh->topology() == meshtags.topology());
const int tdim = mesh->topology()->dim();
const int gdim = mesh->geometry().dim();
const int fdim = tdim - 1;
const int block_size = V->dofmap()->index_map_bs();
std::int32_t size_local = V->dofmap()->index_map->size_local();
const int block_size = V.dofmap()->index_map_bs();
std::int32_t size_local = V.dofmap()->index_map->size_local();

mesh->topology_mutable()->create_connectivity(fdim, tdim);
mesh->topology_mutable()->create_connectivity(tdim, tdim);
Expand Down Expand Up @@ -408,7 +406,7 @@ mpc_data<T> create_contact_slip_condition(
// Note that this function casts the normal array from being potentially
// complex to real valued
std::vector<std::int32_t> dofs(block_size);
std::span<const T> normal_array = nh->x()->array();
std::span<const T> normal_array = nh.x()->array();
const auto largest_normal_component
= [&dofs, block_size, &normal_array, gdim](const std::int32_t block,
std::span<U, 3> normal)
Expand Down Expand Up @@ -449,19 +447,19 @@ mpc_data<T> create_contact_slip_condition(

// Create map from slave dof blocks to a cell containing them
std::vector<std::int32_t> slave_cells = dolfinx_mpc::create_block_to_cell_map(
*mesh->topology(), *V->dofmap(), local_slave_blocks);
*mesh->topology(), *V.dofmap(), local_slave_blocks);
std::vector<U> slave_coordinates;
std::array<std::size_t, 2> coord_shape;
{
std::tie(slave_coordinates, coord_shape)
= dolfinx_mpc::tabulate_dof_coordinates<U>(*V, local_slave_blocks,
= dolfinx_mpc::tabulate_dof_coordinates<U>(V, local_slave_blocks,
slave_cells);
std::vector<std::int32_t> local_cell_collisions
= dolfinx_mpc::find_local_collisions<U>(*mesh, bb_tree,
slave_coordinates, eps2);

auto [basis, basis_shape] = dolfinx_mpc::evaluate_basis_functions<U>(
*V, slave_coordinates, local_cell_collisions);
V, slave_coordinates, local_cell_collisions);
assert(basis_shape.back() == 1);
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
Expand Down Expand Up @@ -582,7 +580,7 @@ mpc_data<T> create_contact_slip_condition(
= dolfinx_mpc::find_local_collisions<U>(*mesh, bb_tree, recv_coords,
eps2);
auto [recv_basis_values, shape] = dolfinx_mpc::evaluate_basis_functions<U>(
*V, recv_coords, remote_cell_collisions);
V, recv_coords, remote_cell_collisions);
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
basis_span(recv_basis_values.data(), shape[0], shape[1]);
Expand Down Expand Up @@ -884,30 +882,30 @@ mpc_data<T> create_contact_slip_condition(
/// collision
template <typename T, std::floating_point U>
mpc_data<T> create_contact_inelastic_condition(
std::shared_ptr<dolfinx::fem::FunctionSpace<U>> V,
const dolfinx::fem::FunctionSpace<U>& V,
dolfinx::mesh::MeshTags<std::int32_t> meshtags, std::int32_t slave_marker,
std::int32_t master_marker, const double eps2 = 1e-20)
std::int32_t master_marker, const U eps2 = 1e-20)
{
dolfinx::common::Timer timer("~MPC: Inelastic condition");

MPI_Comm comm = V->mesh()->comm();
MPI_Comm comm = V.mesh()->comm();
int rank = -1;
MPI_Comm_rank(comm, &rank);

// Extract some const information from function-space
const std::shared_ptr<const dolfinx::common::IndexMap> imap
= V->dofmap()->index_map;
const int tdim = V->mesh()->topology()->dim();
= V.dofmap()->index_map;
const int tdim = V.mesh()->topology()->dim();
const int fdim = tdim - 1;
const int block_size = V->dofmap()->index_map_bs();
std::int32_t size_local = V->dofmap()->index_map->size_local();
const int block_size = V.dofmap()->index_map_bs();
std::int32_t size_local = V.dofmap()->index_map->size_local();

// Create entity permutations needed in evaluate_basis_functions
V->mesh()->topology_mutable()->create_entity_permutations();
V.mesh()->topology_mutable()->create_entity_permutations();
// Create connectivities needed for evaluate_basis_functions and
// select_colliding cells
V->mesh()->topology_mutable()->create_connectivity(fdim, tdim);
V->mesh()->topology_mutable()->create_connectivity(tdim, tdim);
V.mesh()->topology_mutable()->create_connectivity(fdim, tdim);
V.mesh()->topology_mutable()->create_connectivity(tdim, tdim);

std::vector<std::int32_t> slave_blocks
= impl::locate_slave_dofs<U>(V, meshtags, slave_marker);
Expand Down Expand Up @@ -962,19 +960,19 @@ mpc_data<T> create_contact_inelastic_condition(
}

// Create boundingboxtree for master surface
auto facet_to_cell = V->mesh()->topology()->connectivity(fdim, tdim);
auto facet_to_cell = V.mesh()->topology()->connectivity(fdim, tdim);
assert(facet_to_cell);
dolfinx::geometry::BoundingBoxTree<U> bb_tree = impl::create_boundingbox_tree(
*V->mesh(), meshtags, master_marker, std::sqrt(eps2));
*V.mesh(), meshtags, master_marker, std::sqrt(eps2));

// Tabulate slave block coordinates and find colliding cells
std::vector<std::int32_t> slave_cells = dolfinx_mpc::create_block_to_cell_map(
*V->mesh()->topology(), *V->dofmap(), local_blocks);
*V.mesh()->topology(), *V.dofmap(), local_blocks);
std::vector<U> slave_coordinates;
{
std::array<std::size_t, 2> c_shape;
std::tie(slave_coordinates, c_shape)
= dolfinx_mpc::tabulate_dof_coordinates<U>(*V, local_blocks,
= dolfinx_mpc::tabulate_dof_coordinates<U>(V, local_blocks,
slave_cells);
}

Expand All @@ -987,10 +985,10 @@ mpc_data<T> create_contact_inelastic_condition(
std::vector<size_t> collision_to_local;
{
std::vector<std::int32_t> colliding_cells
= dolfinx_mpc::find_local_collisions<U>(*V->mesh(), bb_tree,
= dolfinx_mpc::find_local_collisions<U>(*V.mesh(), bb_tree,
slave_coordinates, eps2);
auto [basis_values, basis_shape] = dolfinx_mpc::evaluate_basis_functions<U>(
*V, slave_coordinates, colliding_cells);
V, slave_coordinates, colliding_cells);
assert(basis_shape.back() == 1);
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
Expand All @@ -1003,7 +1001,7 @@ mpc_data<T> create_contact_inelastic_condition(
{
if (const auto& cell = colliding_cells[i]; cell != -1)
{
auto cell_blocks = V->dofmap()->cell_dofs(cell);
auto cell_blocks = V.dofmap()->cell_dofs(cell);
l_master.reserve(cell_blocks.size());
coeff.reserve(cell_blocks.size());
assert(l_master.empty());
Expand Down Expand Up @@ -1172,10 +1170,10 @@ mpc_data<T> create_contact_inelastic_condition(
collision_block_offsets(indegree);
{
std::vector<std::int32_t> remote_cell_collisions
= dolfinx_mpc::find_local_collisions<U>(*V->mesh(), bb_tree,
recv_coords, eps2);
= dolfinx_mpc::find_local_collisions<U>(*V.mesh(), bb_tree, recv_coords,
eps2);
auto [basis, basis_shape] = dolfinx_mpc::evaluate_basis_functions<U>(
*V, recv_coords, remote_cell_collisions);
V, recv_coords, remote_cell_collisions);
assert(basis_shape.back() == 1);
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
Expand All @@ -1196,7 +1194,7 @@ mpc_data<T> create_contact_inelastic_condition(
= std::vector<std::int32_t>(tdim, 0);
if (const auto& cell = remote_cell_collisions[j]; cell != -1)
{
auto cell_blocks = V->dofmap()->cell_dofs(cell);
auto cell_blocks = V.dofmap()->cell_dofs(cell);
r_master.reserve(cell_blocks.size());
r_coeff.reserve(cell_blocks.size());
assert(r_master.empty());
Expand Down
16 changes: 8 additions & 8 deletions cpp/MultiPointConstraint.h
Original file line number Diff line number Diff line change
Expand Up @@ -34,11 +34,11 @@ class MultiPointConstraint
/// @param[in] offsets Offsets for masters
/// @tparam The floating type of the mesh
MultiPointConstraint(std::shared_ptr<const dolfinx::fem::FunctionSpace<U>> V,
const std::vector<std::int32_t>& slaves,
const std::vector<std::int64_t>& masters,
const std::vector<T>& coeffs,
const std::vector<std::int32_t>& owners,
const std::vector<std::int32_t>& offsets)
std::span<const std::int32_t> slaves,
std::span<const std::int64_t> masters,
std::span<const T> coeffs,
std::span<const std::int32_t> owners,
std::span<const std::int32_t> offsets)
: _slaves(), _is_slave(), _cell_to_slaves_map(), _num_local_slaves(),
_master_map(), _coeff_map(), _owner_map(), _mpc_constants(), _V()
{
Expand All @@ -64,7 +64,7 @@ class MultiPointConstraint
_is_slave = std::move(_slave_data);

// Create a map for cells owned by the process to the slaves
_cell_to_slaves_map = create_cell_to_dofs_map(V, slaves);
_cell_to_slaves_map = create_cell_to_dofs_map(*V, slaves);

// Create adjacency list with all local dofs, where the slave dofs maps to
// its masters
Expand Down Expand Up @@ -117,11 +117,11 @@ class MultiPointConstraint

// Create new function space with extended index map
_V = std::make_shared<const dolfinx::fem::FunctionSpace<U>>(
create_extended_functionspace(V, _master_data, _owner_data));
create_extended_functionspace(*V, _master_data, _owner_data));

// Map global masters to local index in extended function space
std::vector<std::int32_t> masters_local
= map_dofs_global_to_local<U>(_V, _master_data);
= map_dofs_global_to_local<U>(*_V, _master_data);
_master_map = std::make_shared<dolfinx::graph::AdjacencyList<std::int32_t>>(
masters_local, masters_offsets);
}
Expand Down
8 changes: 5 additions & 3 deletions cpp/PeriodicConstraint.h
Original file line number Diff line number Diff line change
Expand Up @@ -129,13 +129,13 @@ dolfinx_mpc::mpc_data<T> _create_periodic_condition(
// Create bounding-box tree over owned cells
std::vector<std::int32_t> r(num_cells_local);
std::iota(r.begin(), r.end(), 0);
dolfinx::geometry::BoundingBoxTree<U> tree(*mesh.get(), tdim, r, 1e-15);
dolfinx::geometry::BoundingBoxTree<U> tree(*mesh.get(), tdim, r, tol);
auto process_tree = tree.create_global_tree(mesh->comm());
auto colliding_bbox_processes
= dolfinx::geometry::compute_collisions<U>(process_tree, mapped_T_b);

std::vector<std::int32_t> local_cell_collisions
= dolfinx_mpc::find_local_collisions<U>(*mesh, tree, mapped_T_b, 1e-20);
= dolfinx_mpc::find_local_collisions<U>(*mesh, tree, mapped_T_b, tol);
dolfinx::common::Timer t0("~~Periodic: Local cell and eval basis");
auto [basis_values, basis_shape] = dolfinx_mpc::evaluate_basis_functions<U>(
V, mapped_T_b, local_cell_collisions);
Expand Down Expand Up @@ -532,6 +532,7 @@ dolfinx_mpc::mpc_data<T> geometrical_condition(
{
std::vector<std::int32_t> slave_blocks
= dolfinx::fem::locate_dofs_geometrical(*V, indicator);

reduced_blocks.reserve(slave_blocks.size());
// Remove blocks in Dirichlet bcs
std::vector<std::int8_t> bc_marker
Expand Down Expand Up @@ -568,7 +569,8 @@ dolfinx_mpc::mpc_data<T> topological_condition(
T scale, bool collapse)
{
std::vector<std::int32_t> entities = meshtag->find(tag);

V->mesh()->topology_mutable()->create_connectivity(
meshtag->dim(), V->mesh()->topology()->dim());
if (collapse)
{
// Locate dofs in sub and parent space
Expand Down
Loading