9#include "dolfinx/mesh/Mesh.h" 
   10#include "dolfinx/mesh/cell_types.h" 
   12#include "dolfinx/refinement/plaza.h" 
   22#include "dolfinx/refinement/option.h" 
   23#include "dolfinx/refinement/utils.h" 
   25namespace dolfinx::refinement::interval
 
   35template <std::
floating_po
int T>
 
   36std::tuple<graph::AdjacencyList<std::int64_t>, std::vector<T>,
 
   37           std::array<std::size_t, 2>, std::optional<std::vector<std::int32_t>>,
 
   38           std::optional<std::vector<std::int8_t>>>
 
   39compute_refinement_data(
const mesh::Mesh<T>& mesh,
 
   40                        std::optional<std::span<const std::int32_t>> cells,
 
   46  if (compute_parent_facet)
 
   47    throw std::runtime_error(
"Parent facet computation not yet supported!");
 
   49  auto topology = mesh.topology();
 
   51  assert(topology->dim() == 1);
 
   52  auto map_c = topology->index_map(1);
 
   58  graph::AdjacencyList<int> cell_ranks = map_c->index_to_dest_ranks();
 
   62  std::vector<int> ranks = cell_ranks.array();
 
   63  std::ranges::sort(ranks);
 
   64  auto to_remove = std::ranges::unique(ranks);
 
   65  ranks.erase(to_remove.begin(), to_remove.end());
 
   68  std::ranges::transform(cell_ranks.array(), cell_ranks.array().begin(),
 
   71                           auto it = std::lower_bound(ranks.begin(),
 
   73                           assert(it != ranks.end() and *it == r);
 
   74                           return std::distance(ranks.begin(), it);
 
   78  std::vector<std::int8_t> refinement_marker(
 
   79      map_c->size_local() + map_c->num_ghosts(), !
cells.has_value());
 
   82  std::vector<std::vector<std::int32_t>> marked_for_update(ranks.size());
 
   85    std::ranges::for_each(
 
   87        [&refinement_marker, &cell_ranks, &marked_for_update](
auto cell)
 
   89          if (!refinement_marker[cell])
 
   91            refinement_marker[cell] = true;
 
   92            for (int rank : cell_ranks.links(cell))
 
   93              marked_for_update[rank].push_back(cell);
 
   99  MPI_Comm neighbor_comm;
 
  100  MPI_Dist_graph_create_adjacent(
 
  101      mesh.comm(), ranks.size(), ranks.data(), MPI_UNWEIGHTED, ranks.size(),
 
  102      ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, 
false, &neighbor_comm);
 
  109                                refinement_marker, *map_c);
 
  113  const auto [new_vertex_map, new_vertex_coords, xshape]
 
  115  MPI_Comm_free(&neighbor_comm);
 
  117  auto c_to_v = mesh.topology()->connectivity(1, 0);
 
  122  const std::int32_t number_of_refined_cells
 
  123      = std::count(refinement_marker.begin(),
 
  124                   std::next(refinement_marker.begin(),
 
  125                             mesh.topology()->index_map(1)->size_local()),
 
  129  std::vector<std::int64_t> global_indices
 
  130      = 
adjust_indices(*mesh.topology()->index_map(0), number_of_refined_cells);
 
  133  std::size_t refined_cell_count
 
  134      = mesh.topology()->index_map(1)->size_local() + number_of_refined_cells;
 
  136  std::vector<std::int64_t> cell_topology;
 
  137  cell_topology.reserve(refined_cell_count * 2);
 
  139  std::optional<std::vector<std::int32_t>> 
parent_cell(std::nullopt);
 
  140  if (compute_parent_cell)
 
  146  for (std::int32_t cell = 0; cell < map_c->size_local(); ++cell)
 
  148    std::span vertices = c_to_v->links(cell);
 
  149    assert(vertices.size() == 2);
 
  153    const std::int64_t a = global_indices[vertices[0]];
 
  154    const std::int64_t b = global_indices[vertices[1]];
 
  155    if (refinement_marker[cell])
 
  159      std::span nv = new_vertex_map.links(cell);
 
  160      assert(nv.size() == 1);
 
  161      const std::int64_t c = nv[0];
 
  164      cell_topology.insert(cell_topology.end(), {a, c, c, b});
 
  166      if (compute_parent_cell)
 
  172      cell_topology.insert(cell_topology.end(), {a, b});
 
  174      if (compute_parent_cell)
 
  179  assert(cell_topology.size() == 2 * refined_cell_count);
 
  180  assert(!compute_parent_cell or 
parent_cell->size() == refined_cell_count);
 
  182  std::vector<std::int32_t> offsets(refined_cell_count + 1);
 
  183  std::ranges::generate(offsets, [i = 0]() 
mutable { 
return 2 * i++; });
 
  185  graph::AdjacencyList cell_adj(std::move(cell_topology), std::move(offsets));
 
  187  return {std::move(cell_adj), std::move(new_vertex_coords), xshape,
 
Functions supporting mesh operations.
 
void cells(la::SparsityPattern &pattern, std::array< std::span< const std::int32_t >, 2 > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:16
 
std::tuple< graph::AdjacencyList< std::int64_t >, std::vector< T >, std::array< std::size_t, 2 > > create_new_vertices(MPI_Comm comm, const graph::AdjacencyList< int > &shared_edges, const mesh::Mesh< T > &mesh, std::span< const std::int8_t > marked_edges)
Add new vertex for each marked edge, and create new_vertex_coordinates and global_edge->new_vertex ma...
Definition utils.h:148
 
constexpr bool option_parent_cell(Option a)
Check if parent_cell flag is set.
Definition option.h:45
 
Option
Options for data to compute during mesh refinement.
Definition option.h:16
 
@ parent_cell
Definition option.h:20
 
constexpr bool option_parent_facet(Option a)
Check if parent_facet flag is set.
Definition option.h:36
 
std::vector< std::int64_t > adjust_indices(const common::IndexMap &map, std::int32_t n)
Given an index map, add "n" extra indices at the end of local range.
Definition utils.cpp:103
 
void update_logical_edgefunction(MPI_Comm comm, const std::vector< std::vector< std::int32_t > > &marked_for_update, std::span< std::int8_t > marked_edges, const common::IndexMap &map)
Communicate edge markers between processes that share edges.
Definition utils.cpp:47