Skip to content

Commit

Permalink
Fixes to be compatible with DOLFINx main (april 20, 2023) (#60)
Browse files Browse the repository at this point in the history
* Fix integration domain handling (FEniCS/dolfinx#2617)

* Fix lifting
  • Loading branch information
jorgensd committed Apr 20, 2023
1 parent 2af8c96 commit 653b5c3
Show file tree
Hide file tree
Showing 5 changed files with 19 additions and 13 deletions.
10 changes: 6 additions & 4 deletions cpp/assemble_matrix.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -263,7 +263,7 @@ void assemble_exterior_facets(
const std::span<const std::int32_t>&,
const std::span<const T>)>& mat_add_values,
const dolfinx::mesh::Mesh<double>& mesh,
const std::span<const std::int32_t>& facets,
std::span<const std::int32_t> facets,
const std::function<void(const std::span<T>&,
const std::span<const std::uint32_t>&,
std::int32_t, int)>& apply_dof_transformation,
Expand Down Expand Up @@ -396,7 +396,7 @@ void assemble_cells_impl(
const std::span<const std::int32_t>&,
const std::span<const T>)>& mat_add_values,
const dolfinx::mesh::Geometry<double>& geometry,
const std::vector<std::int32_t>& active_cells,
std::span<const std::int32_t> active_cells,
std::function<void(std::span<T>, const std::span<const std::uint32_t>,
const std::int32_t, const int)>
apply_dof_transformation,
Expand Down Expand Up @@ -571,7 +571,8 @@ void assemble_matrix_impl(
const auto& fn = a.kernel(dolfinx::fem::IntegralType::cell, i);
const auto& [coeffs, cstride]
= coefficients.at({dolfinx::fem::IntegralType::cell, i});
const std::vector<std::int32_t>& active_cells = a.cell_domains(i);
std::span<const std::int32_t> active_cells
= a.domain(dolfinx::fem::IntegralType::cell, i);
assemble_cells_impl<T>(
mat_add_block_values, mat_add_values, mesh->geometry(), active_cells,
apply_dof_transformation, dofs0, bs0,
Expand All @@ -584,7 +585,8 @@ void assemble_matrix_impl(
const auto& fn = a.kernel(dolfinx::fem::IntegralType::exterior_facet, i);
const auto& [coeffs, cstride]
= coefficients.at({dolfinx::fem::IntegralType::exterior_facet, i});
const std::vector<std::int32_t>& facets = a.exterior_facet_domains(i);
std::span<const std::int32_t> facets
= a.domain(dolfinx::fem::IntegralType::exterior_facet, i);
assemble_exterior_facets<T>(mat_add_block_values, mat_add_values, *mesh,
facets, apply_dof_transformation, dofs0, bs0,
apply_dof_transformation_to_transpose, dofs1,
Expand Down
7 changes: 4 additions & 3 deletions cpp/assemble_vector.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -167,7 +167,8 @@ void _assemble_vector(
};

// Assemble over all active cells
const std::vector<std::int32_t>& active_cells = L.cell_domains(i);
std::span<const std::int32_t> active_cells
= L.domain(dolfinx::fem::IntegralType::cell, i);
_assemble_entities_impl<T, 1>(b, active_cells, dofs, bs, mpc, fetch_cell,
assemble_local_cell_vector);
}
Expand Down Expand Up @@ -214,8 +215,8 @@ void _assemble_vector(
};

// Assemble over all active cells
const std::vector<std::int32_t>& active_facets
= L.exterior_facet_domains(i);
std::span<const std::int32_t> active_facets
= L.domain(dolfinx::fem::IntegralType::exterior_facet, i);
_assemble_entities_impl<T, 2>(b, active_facets, dofs, bs, mpc, fetch_cell,
assemble_local_exterior_facet_vector);
}
Expand Down
7 changes: 4 additions & 3 deletions cpp/lifting.h
Original file line number Diff line number Diff line change
Expand Up @@ -278,7 +278,8 @@ void apply_lifting(
}
};
// Assemble over all active cells
const std::vector<std::int32_t>& cells = a->cell_domains(i);
std::span<const std::int32_t> cells
= a->domain(dolfinx::fem::IntegralType::cell, i);
lift_bc_entities<T, 1>(b, cells, dofmap0, dofmap1, bs0, bs1, bc_values1,
bc_markers1, mpc1, fetch_cells, lift_bcs_cell);
}
Expand Down Expand Up @@ -363,8 +364,8 @@ void apply_lifting(
};

// Assemble over all active cells
const std::vector<std::int32_t>& active_facets
= a->exterior_facet_domains(i);
std::span<const std::int32_t> active_facets
= a->domain(dolfinx::fem::IntegralType::exterior_facet, i);
impl::lift_bc_entities<T, 2>(b, active_facets, dofmap0, dofmap1, bs0, bs1,
bc_values1, bc_markers1, mpc1, fetch_cell,
lift_bc_exterior_facet);
Expand Down
1 change: 1 addition & 0 deletions cpp/utils.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,7 @@
#include <algorithm>
#include <basix/mdspan.hpp>
#include <dolfinx/geometry/utils.h>
#include <dolfinx/la/petsc.h>
#include <dolfinx/mesh/Mesh.h>
#include <dolfinx/mesh/utils.h>

Expand Down
7 changes: 4 additions & 3 deletions cpp/utils.h
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@
#include <dolfinx/graph/AdjacencyList.h>
#include <dolfinx/la/SparsityPattern.h>
#include <dolfinx/la/petsc.h>
#include <dolfinx/mesh/MeshTags.h>
#include <span>

namespace impl
Expand Down Expand Up @@ -285,15 +286,15 @@ void build_standard_pattern(dolfinx::la::SparsityPattern& pattern,
case dolfinx::fem::IntegralType::cell:
for (int id : ids)
{
const std::vector<std::int32_t>& cells = a.cell_domains(id);
std::span<const std::int32_t> cells = a.domain(type, id);
dolfinx::fem::sparsitybuild::cells(pattern, cells,
{{dofmaps[0], dofmaps[1]}});
}
break;
case dolfinx::fem::IntegralType::interior_facet:
for (int id : ids)
{
const std::vector<std::int32_t>& facets = a.interior_facet_domains(id);
std::span<const std::int32_t> facets = a.domain(type, id);
std::vector<std::int32_t> f;
f.reserve(facets.size() / 2);
for (std::size_t i = 0; i < facets.size(); i += 4)
Expand All @@ -305,7 +306,7 @@ void build_standard_pattern(dolfinx::la::SparsityPattern& pattern,
case dolfinx::fem::IntegralType::exterior_facet:
for (int id : ids)
{
const std::vector<std::int32_t>& facets = a.exterior_facet_domains(id);
std::span<const std::int32_t> facets = a.domain(type, id);
std::vector<std::int32_t> cells;
cells.reserve(facets.size() / 2);
for (std::size_t i = 0; i < facets.size(); i += 2)
Expand Down

0 comments on commit 653b5c3

Please sign in to comment.