# include <deal.II/fe/mapping.h>
# include <deal.II/fe/mapping_q1.h>
-# include <deal.II/lac/constraint_matrix.h>
+# include <deal.II/lac/affine_constraints.h>
# include <deal.II/lac/cuda_vector.h>
# include <cuda_runtime_api.h>
/**
* Extracts the information needed to perform loops over cells. The
- * DoFHandler and ConstraintMatrix describe the layout of degrees of
- * freedom, the DoFHandler and the mapping describe the transformation from
- * unit to real cell, and the finite element underlying the DoFHandler
- * together with the quadrature formula describe the local operations.
+ * DoFHandler and AffineConstraints objects describe the layout of
+ * degrees of freedom, the DoFHandler and the mapping describe the
+ * transformation from unit to real cell, and the finite element
+ * underlying the DoFHandler together with the quadrature formula
+ * describe the local operations.
*/
void
- reinit(const Mapping<dim> & mapping,
- const DoFHandler<dim> & dof_handler,
- const ConstraintMatrix &constraints,
- const Quadrature<1> & quad,
- const AdditionalData additional_data = AdditionalData());
+ reinit(const Mapping<dim> & mapping,
+ const DoFHandler<dim> & dof_handler,
+ const AffineConstraints<Number> &constraints,
+ const Quadrature<1> & quad,
+ const AdditionalData additional_data = AdditionalData());
/**
* Initializes the data structures. Same as above but using a Q1 mapping.
*/
void
- reinit(const DoFHandler<dim> & dof_handler,
- const ConstraintMatrix &constraints,
- const Quadrature<1> & quad,
- const AdditionalData AdditionalData = AdditionalData());
+ reinit(const DoFHandler<dim> & dof_handler,
+ const AffineConstraints<Number> &constraints,
+ const Quadrature<1> & quad,
+ const AdditionalData AdditionalData = AdditionalData());
/**
* Return the Data structure associated with @p color.
- template <int dim>
+ template <int dim, typename number>
std::vector<types::global_dof_index>
get_conflict_indices(
const FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
- & cell,
- const ConstraintMatrix &constraints)
+ & cell,
+ const AffineConstraints<number> &constraints)
{
std::vector<types::global_dof_index> local_dof_indices(
cell->get_fe().dofs_per_cell);
template <int dim, typename Number>
void
- MatrixFree<dim, Number>::reinit(const Mapping<dim> & mapping,
- const DoFHandler<dim> & dof_handler,
- const ConstraintMatrix &constraints,
- const Quadrature<1> & quad,
- const AdditionalData additional_data)
+ MatrixFree<dim, Number>::reinit(const Mapping<dim> & mapping,
+ const DoFHandler<dim> & dof_handler,
+ const AffineConstraints<Number> &constraints,
+ const Quadrature<1> & quad,
+ const AdditionalData additional_data)
{
if (typeid(Number) == typeid(double))
cudaDeviceSetSharedMemConfig(cudaSharedMemBankSizeEightByte);
CellFilter end(IteratorFilters::LocallyOwnedCell(), dof_handler.end());
typedef std::function<std::vector<types::global_dof_index>(
CellFilter const &)>
- fun_type;
- const fun_type &fun =
- static_cast<fun_type>(std::bind(&internal::get_conflict_indices<dim>,
- std::placeholders::_1,
- constraints));
+ fun_type;
+ const auto fun = [&](const CellFilter &filter) {
+ return internal::get_conflict_indices<dim, Number>(filter, constraints);
+ };
std::vector<std::vector<CellFilter>> graph =
GraphColoring::make_graph_coloring(begin, end, fun);