Skip to content

Commit

Permalink
Allow AffineConstraints argument to MGTools::make_sparsity_pattern
Browse files Browse the repository at this point in the history
  • Loading branch information
kronbichler committed Aug 11, 2020
1 parent 53bf694 commit f5bc774
Show file tree
Hide file tree
Showing 3 changed files with 42 additions and 15 deletions.
23 changes: 16 additions & 7 deletions include/deal.II/multigrid/mg_tools.h
Original file line number Diff line number Diff line change
Expand Up @@ -75,16 +75,25 @@ namespace MGTools
* Write the sparsity structure of the matrix belonging to the specified @p
* level. The sparsity pattern is not compressed, so before creating the
* actual matrix you have to compress the matrix yourself, using
* <tt>SparseMatrixStruct::compress()</tt>.
* <tt>SparsityPatternType::compress()</tt>.
*
* There is no need to consider hanging nodes here, since only one level is
* considered.
* The optional AffineConstraints argument allows to define constraints of
* the level matrices like Dirichlet boundary conditions. Note that there is
* need to consider hanging nodes on the typical level matrices, since only
* one level is considered. See DoFTools::make_sparsity_pattern() for more
* details about the arguments.
*/
template <int dim, int spacedim, typename SparsityPatternType>
template <int dim,
int spacedim,
typename SparsityPatternType,
typename number = double>
void
make_sparsity_pattern(const DoFHandler<dim, spacedim> &dof_handler,
SparsityPatternType & sparsity,
const unsigned int level);
make_sparsity_pattern(
const DoFHandler<dim, spacedim> &dof_handler,
SparsityPatternType & sparsity,
const unsigned int level,
const AffineConstraints<number> &constraints = AffineConstraints<number>(),
const bool keep_constrained_dofs = true);

/**
* Make a sparsity pattern including fluxes of discontinuous Galerkin
Expand Down
16 changes: 10 additions & 6 deletions source/multigrid/mg_tools.cc
Original file line number Diff line number Diff line change
Expand Up @@ -567,11 +567,16 @@ namespace MGTools



template <int dim, int spacedim, typename SparsityPatternType>
template <int dim,
int spacedim,
typename SparsityPatternType,
typename number>
void
make_sparsity_pattern(const DoFHandler<dim, spacedim> &dof,
SparsityPatternType & sparsity,
const unsigned int level)
const unsigned int level,
const AffineConstraints<number> &constraints,
const bool keep_constrained_dofs)
{
const types::global_dof_index n_dofs = dof.n_dofs(level);
(void)n_dofs;
Expand All @@ -592,10 +597,9 @@ namespace MGTools
dof.get_triangulation().locally_owned_subdomain())
{
cell->get_mg_dof_indices(dofs_on_this_cell);
// make sparsity pattern for this cell
for (unsigned int i = 0; i < dofs_per_cell; ++i)
for (unsigned int j = 0; j < dofs_per_cell; ++j)
sparsity.add(dofs_on_this_cell[i], dofs_on_this_cell[j]);
constraints.add_entries_local_to_global(dofs_on_this_cell,
sparsity,
keep_constrained_dofs);
}
}

Expand Down
18 changes: 16 additions & 2 deletions source/multigrid/mg_tools.inst.in
Original file line number Diff line number Diff line change
Expand Up @@ -25,10 +25,24 @@ for (PATTERN : SPARSITY_PATTERNS; deal_II_dimension : DIMENSIONS;
template void
make_sparsity_pattern<deal_II_dimension,
deal_II_space_dimension,
PATTERN>(
PATTERN,
double>(
const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
PATTERN &,
const unsigned int);
const unsigned int,
const AffineConstraints<double> &,
const bool);

template void
make_sparsity_pattern<deal_II_dimension,
deal_II_space_dimension,
PATTERN,
float>(
const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
PATTERN &,
const unsigned int,
const AffineConstraints<float> &,
const bool);
#endif
\}
}
Expand Down

0 comments on commit f5bc774

Please sign in to comment.