--- /dev/null
+Deprecated: All `MatrixFree::reinit()` functions without `Mapping` as argument
+have been deprecated. Use the functions that explicit provide the Mapping instead.
+<br>
+(Peter Munch, 2020/11/11)
FE_Q<dim> fe;
DoFHandler<dim> dof_handler;
+ MappingQ1<dim> mapping;
+
AffineConstraints<double> constraints;
using SystemMatrixType =
LaplaceOperator<dim, degree_finite_element, double>;
constraints.clear();
constraints.reinit(locally_relevant_dofs);
DoFTools::make_hanging_node_constraints(dof_handler, constraints);
- VectorTools::interpolate_boundary_values(dof_handler,
- 0,
- Functions::ZeroFunction<dim>(),
- constraints);
+ VectorTools::interpolate_boundary_values(
+ mapping, dof_handler, 0, Functions::ZeroFunction<dim>(), constraints);
constraints.close();
setup_time += time.wall_time();
time_details << "Distribute DoFs & B.C. (CPU/wall) " << time.cpu_time()
(update_gradients | update_JxW_values | update_quadrature_points);
std::shared_ptr<MatrixFree<dim, double>> system_mf_storage(
new MatrixFree<dim, double>());
- system_mf_storage->reinit(dof_handler,
+ system_mf_storage->reinit(mapping,
+ dof_handler,
constraints,
QGauss<1>(fe.degree + 1),
additional_data);
additional_data.mg_level = level;
std::shared_ptr<MatrixFree<dim, float>> mg_mf_storage_level(
new MatrixFree<dim, float>());
- mg_mf_storage_level->reinit(dof_handler,
+ mg_mf_storage_level->reinit(mapping,
+ dof_handler,
level_constraints,
QGauss<1>(fe.degree + 1),
additional_data);
solution.update_ghost_values();
data_out.attach_dof_handler(dof_handler);
data_out.add_data_vector(solution, "solution");
- data_out.build_patches();
+ data_out.build_patches(mapping);
DataOutBase::VtkFlags flags;
flags.compression_level = DataOutBase::VtkFlags::best_speed;
#else
Triangulation<dim> triangulation;
#endif
- FE_Q<dim> fe;
- DoFHandler<dim> dof_handler;
+ FE_Q<dim> fe;
+ DoFHandler<dim> dof_handler;
+
+ MappingQ1<dim> mapping;
+
AffineConstraints<double> constraints;
IndexSet locally_relevant_dofs;
additional_data.tasks_parallel_scheme =
MatrixFree<dim>::AdditionalData::TasksParallelScheme::partition_partition;
- matrix_free_data.reinit(dof_handler,
+ matrix_free_data.reinit(mapping,
+ dof_handler,
constraints,
QGaussLobatto<1>(fe_degree + 1),
additional_data);
Vector<float> norm_per_cell(triangulation.n_active_cells());
solution.update_ghost_values();
- VectorTools::integrate_difference(dof_handler,
+ VectorTools::integrate_difference(mapping,
+ dof_handler,
solution,
Functions::ZeroFunction<dim>(),
norm_per_cell,
data_out.attach_dof_handler(dof_handler);
data_out.add_data_vector(solution, "solution");
- data_out.build_patches();
+ data_out.build_patches(mapping);
data_out.write_vtu_with_pvtu_record(
"./", "solution", timestep_number, MPI_COMM_WORLD, 3);
// get later consumed by the SineGordonOperation::apply() function. Next,
// an instance of the <code> SineGordonOperation class </code> based on
// the finite element degree specified at the top of this file is set up.
- VectorTools::interpolate(dof_handler,
+ VectorTools::interpolate(mapping,
+ dof_handler,
InitialCondition<dim>(1, time),
solution);
- VectorTools::interpolate(dof_handler,
+ VectorTools::interpolate(mapping,
+ dof_handler,
InitialCondition<dim>(1, time - time_step),
old_solution);
output_results(0);
(update_gradients | update_JxW_values | update_quadrature_points);
std::shared_ptr<MatrixFree<dim, double>> mf_storage =
std::make_shared<MatrixFree<dim, double>>();
- mf_storage->reinit(dof_handler,
+ mf_storage->reinit(mapping,
+ dof_handler,
constraints,
QGauss<1>(degree + 1),
additional_data);
additional_data.mg_level = level;
std::shared_ptr<MatrixFree<dim, float>> mf_storage_level(
new MatrixFree<dim, float>());
- mf_storage_level->reinit(dof_handler,
+ mf_storage_level->reinit(mapping,
+ dof_handler,
level_constraints,
QGauss<1>(degree + 1),
additional_data);
FE_DGQHermite<dim> fe;
DoFHandler<dim> dof_handler;
+ MappingQ1<dim> mapping;
+
using SystemMatrixType = LaplaceOperator<dim, fe_degree, double>;
SystemMatrixType system_matrix;
update_quadrature_points);
const auto system_mf_storage =
std::make_shared<MatrixFree<dim, double>>();
- system_mf_storage->reinit(dof_handler,
- dummy,
- QGauss<1>(fe.degree + 1),
- additional_data);
+ system_mf_storage->reinit(
+ mapping, dof_handler, dummy, QGauss<1>(fe.degree + 1), additional_data);
system_matrix.initialize(system_mf_storage);
}
additional_data.mg_level = level;
const auto mg_mf_storage_level =
std::make_shared<MatrixFree<dim, float>>();
- mg_mf_storage_level->reinit(dof_handler,
+ mg_mf_storage_level->reinit(mapping,
+ dof_handler,
dummy,
QGauss<1>(fe.degree + 1),
additional_data);
void LaplaceProblem<dim, fe_degree>::analyze_results() const
{
Vector<float> error_per_cell(triangulation.n_active_cells());
- VectorTools::integrate_difference(dof_handler,
+ VectorTools::integrate_difference(mapping,
+ dof_handler,
solution,
Solution<dim>(),
error_per_cell,
/**
* Initializes the data structures. Same as above, but using a $Q_1$
* mapping.
+ *
+ * @deprecated Use the overload taking a Mapping object instead.
*/
template <typename QuadratureType, typename number2>
- void
+ DEAL_II_DEPRECATED void
reinit(const DoFHandler<dim> & dof_handler,
const AffineConstraints<number2> &constraint,
const QuadratureType & quad,
/**
* Initializes the data structures. Same as above, but using a $Q_1$
* mapping.
+ *
+ * @deprecated Use the overload taking a Mapping object instead.
*/
template <typename QuadratureType, typename number2>
- void
+ DEAL_II_DEPRECATED void
reinit(const std::vector<const DoFHandler<dim> *> & dof_handler,
const std::vector<const AffineConstraints<number2> *> &constraint,
const std::vector<QuadratureType> & quad,
* @deprecated Use the overload taking a DoFHandler object instead.
*/
template <typename QuadratureType, typename number2>
- void
+ DEAL_II_DEPRECATED void
reinit(const std::vector<const hp::DoFHandler<dim> *> & dof_handler_hp,
const std::vector<const AffineConstraints<number2> *> &constraint,
const std::vector<QuadratureType> & quad,
/**
* Initializes the data structures. Same as above, but using a $Q_1$
* mapping.
+ *
+ * @deprecated Use the overload taking a Mapping object instead.
*/
template <typename QuadratureType, typename number2>
- void
+ DEAL_II_DEPRECATED void
reinit(const std::vector<const DoFHandler<dim> *> & dof_handler,
const std::vector<const AffineConstraints<number2> *> &constraint,
const QuadratureType & quad,