1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62
|
// Copyright (C) 2007-2023 Garth N. Wells
//
// This file is part of DOLFINx (https://www.fenicsproject.org)
//
// SPDX-License-Identifier: LGPL-3.0-or-later
#pragma once
#include <array>
#include <cstdint>
#include <functional>
#include <span>
namespace dolfinx::la
{
class SparsityPattern;
}
namespace dolfinx::fem
{
class DofMap;
/// Support for building sparsity patterns from degree-of-freedom maps.
namespace sparsitybuild
{
/// @brief Iterate over cells and insert entries into sparsity pattern.
///
/// Inserts the rectangular blocks of indices `dofmap[0][cells[0][i]] x
/// dofmap[1][cells[1][i]]` into the sparsity pattern, i.e. entries
/// `(dofmap[0][cells[0][i]][k0], dofmap[0][cells[0][i]][k1])` will
/// appear in the sparsity pattern.
///
/// @param pattern Sparsity pattern to insert into.
/// @param cells Lists of cells to iterate over. `cells[0]` and
/// `cells[1]` must have the same size.
/// @param dofmaps Dofmaps to used in building the sparsity pattern.
/// @note The sparsity pattern is not finalised.
void cells(la::SparsityPattern& pattern,
std::array<std::span<const std::int32_t>, 2> cells,
std::array<std::reference_wrapper<const DofMap>, 2> dofmaps);
/// @brief Iterate over interior facets and insert entries into sparsity
/// pattern.
///
/// Inserts the rectangular blocks of indices `[dofmap[0][cell0],
/// dofmap[0][cell1]] x [dofmap[1][cell0] + dofmap[1][cell1]]` where
/// `cell0` and `cell1` are the two cells attached to a facet.
///
/// @param[in,out] pattern Sparsity pattern to insert into.
/// @param[in] cells Cells to index into each dofmap. `cells[i]` is a
/// list of `(cell0, cell1)` pairs for each interior facet to index into
/// `dofmap[i]`. `cells[0]` and `cells[1]` must have the same size.
/// @param[in] dofmaps Dofmaps to use in building the sparsity pattern.
///
/// @note The sparsity pattern is not finalised.
void interior_facets(
la::SparsityPattern& pattern,
std::array<std::span<const std::int32_t>, 2> cells,
std::array<std::reference_wrapper<const DofMap>, 2> dofmaps);
} // namespace sparsitybuild
} // namespace dolfinx::fem
|