Skip to content

Commit

Permalink
Merge pull request #15906 from bangerth/default-2
Browse files Browse the repository at this point in the history
Add a default constructor to AffineConstraints. (And be more terse about default arguments.)
  • Loading branch information
kronbichler committed Aug 23, 2023
2 parents 4ddada2 + e3dc646 commit dfb4350
Show file tree
Hide file tree
Showing 4 changed files with 52 additions and 29 deletions.
6 changes: 3 additions & 3 deletions include/deal.II/dofs/dof_tools.h
Original file line number Diff line number Diff line change
Expand Up @@ -424,7 +424,7 @@ namespace DoFTools
make_sparsity_pattern(
const DoFHandler<dim, spacedim> &dof_handler,
SparsityPatternBase & sparsity_pattern,
const AffineConstraints<number> &constraints = AffineConstraints<number>(),
const AffineConstraints<number> &constraints = {},
const bool keep_constrained_dofs = true,
const types::subdomain_id subdomain_id = numbers::invalid_subdomain_id);

Expand Down Expand Up @@ -499,7 +499,7 @@ namespace DoFTools
const DoFHandler<dim, spacedim> &dof_handler,
const Table<2, Coupling> & coupling,
SparsityPatternBase & sparsity_pattern,
const AffineConstraints<number> &constraints = AffineConstraints<number>(),
const AffineConstraints<number> &constraints = {},
const bool keep_constrained_dofs = true,
const types::subdomain_id subdomain_id = numbers::invalid_subdomain_id);

Expand Down Expand Up @@ -1378,7 +1378,7 @@ namespace DoFTools
const std::function<
bool(const typename DoFHandler<dim, spacedim>::active_cell_iterator &)>
& predicate,
const AffineConstraints<number> &constraints = AffineConstraints<number>());
const AffineConstraints<number> &constraints = {});

/**
* Extract a vector that represents the constant modes of the DoFHandler for
Expand Down
41 changes: 33 additions & 8 deletions include/deal.II/lac/affine_constraints.h
Original file line number Diff line number Diff line change
Expand Up @@ -547,24 +547,38 @@ class AffineConstraints : public Subscriptor
right_object_wins
};

/**
* Default constructor.
*
* This constructor sets up an object that stores all constraints eventually
* given to it. (This is in contrast to the following constructor that
* restricts which degrees of freedom we should care about -- for example,
* in a parallel context one would typically only store constraints
* for the locally relevant DoFs, see
* @ref GlossLocallyRelevantDof .)
* Consequently, this constructor creates internal data structures for
* <i>all</i> possible indices, leading to memory consumption on every
* processor that is proportional to the <i>overall</i> size of the
* problem, not just proportional to the size of the portion of the
* overall problem that is handled by the current processor. Calling
* this constructor in a parallel program is a common bug: Call the
* other one, with an index set, instead.
*/
AffineConstraints();

/**
* Constructor. The supplied IndexSet defines which indices might be
* constrained inside this AffineConstraints container. In a calculation
* with a DoFHandler object based on parallel::distributed::Triangulation
* or parallel::shared::Triangulation, one should use the set of locally
* relevant dofs (see
* relevant DoFs (see
* @ref GlossLocallyRelevantDof).
*
* The given IndexSet allows the AffineConstraints container to save
* memory by just not caring about degrees of freedom that are not of
* importance to the current processor. Alternatively, if no such
* IndexSet is provided, internal data structures for <i>all</i> possible
* indices will be created, leading to memory consumption on every
* processor that is proportional to the <i>overall</i> size of the
* problem, not just proportional to the size of the portion of the
* overall problem that is handled by the current processor.
* importance to the current processor.
*/
explicit AffineConstraints(const IndexSet &local_constraints = IndexSet());
explicit AffineConstraints(const IndexSet &local_constraints);

/**
* Copy constructor
Expand Down Expand Up @@ -2022,6 +2036,13 @@ class AffineConstraints : public Subscriptor

/* ---------------- template and inline functions ----------------- */

template <typename number>
inline AffineConstraints<number>::AffineConstraints()
: AffineConstraints<number>(IndexSet())
{}



template <typename number>
inline AffineConstraints<number>::AffineConstraints(
const IndexSet &local_constraints)
Expand All @@ -2035,6 +2056,8 @@ inline AffineConstraints<number>::AffineConstraints(
local_lines.compress();
}



template <typename number>
inline AffineConstraints<number>::AffineConstraints(
const AffineConstraints &affine_constraints)
Expand All @@ -2045,6 +2068,8 @@ inline AffineConstraints<number>::AffineConstraints(
, sorted(affine_constraints.sorted)
{}



template <typename number>
inline void
AffineConstraints<number>::add_line(const size_type line_n)
Expand Down
22 changes: 10 additions & 12 deletions include/deal.II/multigrid/mg_tools.h
Original file line number Diff line number Diff line change
Expand Up @@ -88,12 +88,11 @@ namespace MGTools
*/
template <int dim, int spacedim, typename number = double>
void
make_sparsity_pattern(
const DoFHandler<dim, spacedim> &dof_handler,
SparsityPatternBase & sparsity,
const unsigned int level,
const AffineConstraints<number> &constraints = AffineConstraints<number>(),
const bool keep_constrained_dofs = true);
make_sparsity_pattern(const DoFHandler<dim, spacedim> &dof_handler,
SparsityPatternBase & sparsity,
const unsigned int level,
const AffineConstraints<number> &constraints = {},
const bool keep_constrained_dofs = true);

/**
* Make a sparsity pattern including fluxes of discontinuous Galerkin
Expand All @@ -105,12 +104,11 @@ namespace MGTools
*/
template <int dim, int spacedim, typename number = double>
void
make_flux_sparsity_pattern(
const DoFHandler<dim, spacedim> &dof_handler,
SparsityPatternBase & sparsity,
const unsigned int level,
const AffineConstraints<number> &constraints = AffineConstraints<number>(),
const bool keep_constrained_dofs = true);
make_flux_sparsity_pattern(const DoFHandler<dim, spacedim> &dof_handler,
SparsityPatternBase & sparsity,
const unsigned int level,
const AffineConstraints<number> &constraints = {},
const bool keep_constrained_dofs = true);


/**
Expand Down
12 changes: 6 additions & 6 deletions include/deal.II/non_matching/coupling.h
Original file line number Diff line number Diff line change
Expand Up @@ -97,8 +97,8 @@ namespace NonMatching
const DoFHandler<dim1, spacedim> &immersed_dh,
const Quadrature<dim1> & quad,
SparsityPatternBase & sparsity,
const AffineConstraints<number> & constraints = AffineConstraints<number>(),
const ComponentMask & space_comps = {},
const AffineConstraints<number> & constraints = {},
const ComponentMask & space_comps = {},
const ComponentMask & immersed_comps = {},
const Mapping<dim0, spacedim> & space_mapping =
StaticMappingQ1<dim0, spacedim>::mapping,
Expand All @@ -121,10 +121,10 @@ namespace NonMatching
const DoFHandler<dim1, spacedim> & immersed_dh,
const Quadrature<dim1> & quad,
SparsityPatternBase & sparsity,
const AffineConstraints<number> &constraints = AffineConstraints<number>(),
const ComponentMask & space_comps = {},
const ComponentMask & immersed_comps = {},
const Mapping<dim1, spacedim> & immersed_mapping =
const AffineConstraints<number> & constraints = {},
const ComponentMask & space_comps = {},
const ComponentMask & immersed_comps = {},
const Mapping<dim1, spacedim> & immersed_mapping =
StaticMappingQ1<dim1, spacedim>::mapping,
const AffineConstraints<number> &immersed_constraints =
AffineConstraints<number>());
Expand Down

0 comments on commit dfb4350

Please sign in to comment.