const DoFHandler<dim> &dof_handler,
const MGConstrainedDoFs &mg_constrained_dofs,
const typename FunctionMap<dim>::type &dirichlet_boundary,
- const unsigned int level = numbers::invalid_unsigned_int)
+ const unsigned int level,
+ const bool threaded)
{
const QGauss<1> quad (n_q_points_1d);
typename MatrixFree<dim,number>::AdditionalData addit_data;
- addit_data.tasks_parallel_scheme = MatrixFree<dim,number>::AdditionalData::none;
+ if (threaded)
+ addit_data.tasks_parallel_scheme = MatrixFree<dim,number>::AdditionalData::partition_partition;
+ else
+ addit_data.tasks_parallel_scheme = MatrixFree<dim,number>::AdditionalData::none;
addit_data.tasks_block_size = 3;
addit_data.level_mg_handler = level;
addit_data.mpi_communicator = MPI_COMM_WORLD;
template <int dim, int fe_degree, int n_q_points_1d, typename number>
-void do_test (const DoFHandler<dim> &dof)
+void do_test (const DoFHandler<dim> &dof, const bool threaded)
{
deallog << "Testing " << dof.get_fe().get_name();
deallog << std::endl;
MappingQ<dim> mapping(fe_degree+1);
LaplaceOperator<dim,fe_degree,n_q_points_1d,double> fine_matrix;
fine_matrix.initialize(mapping, dof, mg_constrained_dofs, dirichlet_boundary,
- numbers::invalid_unsigned_int);
+ numbers::invalid_unsigned_int, threaded);
LinearAlgebra::distributed::Vector<double> in, sol;
fine_matrix.initialize_dof_vector(in);
for (unsigned int level = 0; level<dof.get_triangulation().n_global_levels(); ++level)
{
mg_matrices[level].initialize(mapping, dof, mg_constrained_dofs,
- dirichlet_boundary, level);
+ dirichlet_boundary, level, threaded);
}
MGLevelObject<MGInterfaceMatrix<LevelMatrixType> > mg_interface_matrices;
mg_interface_matrices.resize(0, dof.get_triangulation().n_global_levels()-1);
{
// avoid output from inner (coarse-level) solver
- deallog.depth_file(2);
+ deallog.depth_file(3);
ReductionControl control(30, 1e-20, 1e-7);
SolverCG<LinearAlgebra::distributed::Vector<double> > solver(control);
solver.solve(fine_matrix, sol, in, preconditioner);
for (typename Triangulation<dim>::active_cell_iterator cell=tria.begin_active();
cell != tria.end(); ++cell)
if (cell->is_locally_owned() &&
- (cell->center().norm() < 0.5 && (cell->level() < 5 ||
- cell->center().norm() > 0.45)
+ ((cell->center().norm() < 0.5 && (cell->level() < 5 ||
+ cell->center().norm() > 0.45))
||
(dim == 2 && cell->center().norm() > 1.2)))
cell->set_refine_flag();
dof.distribute_dofs(fe);
dof.distribute_mg_dofs(fe);
- do_test<dim, fe_degree, fe_degree+1, Number> (dof);
+ deallog.push("nothread");
+ do_test<dim, fe_degree, fe_degree+1, Number> (dof, false);
+ deallog.pop();
+ deallog.push("threaded");
+ do_test<dim, fe_degree, fe_degree+1, Number> (dof, true);
+ deallog.pop();
}
}
-DEAL::Testing FE_Q<2>(1)
-DEAL::Number of degrees of freedom: 507
-DEAL:cg::Starting value 21.9317
-DEAL:cg::Convergence step 5 value 7.81205e-07
-DEAL::Testing FE_Q<2>(1)
-DEAL::Number of degrees of freedom: 881
-DEAL:cg::Starting value 27.7669
-DEAL:cg::Convergence step 6 value 9.24482e-07
-DEAL::Testing FE_Q<2>(1)
-DEAL::Number of degrees of freedom: 2147
-DEAL:cg::Starting value 43.2551
-DEAL:cg::Convergence step 6 value 1.14520e-06
-DEAL::Testing FE_Q<2>(1)
-DEAL::Number of degrees of freedom: 6517
-DEAL:cg::Starting value 76.9220
-DEAL:cg::Convergence step 6 value 1.55925e-06
-DEAL::Testing FE_Q<2>(3)
-DEAL::Number of degrees of freedom: 4259
-DEAL:cg::Starting value 64.2573
-DEAL:cg::Convergence step 5 value 1.63117e-06
-DEAL::Testing FE_Q<2>(3)
-DEAL::Number of degrees of freedom: 7457
-DEAL:cg::Starting value 83.1084
-DEAL:cg::Convergence step 5 value 7.69301e-06
-DEAL::Testing FE_Q<2>(3)
-DEAL::Number of degrees of freedom: 18535
-DEAL:cg::Starting value 130.977
-DEAL:cg::Convergence step 5 value 6.04333e-06
-DEAL::Testing FE_Q<3>(1)
-DEAL::Number of degrees of freedom: 186
-DEAL:cg::Starting value 12.3693
-DEAL:cg::Convergence step 3 value 8.40297e-08
-DEAL::Testing FE_Q<3>(1)
-DEAL::Number of degrees of freedom: 648
-DEAL:cg::Starting value 21.1424
-DEAL:cg::Convergence step 6 value 1.05183e-07
-DEAL::Testing FE_Q<3>(1)
-DEAL::Number of degrees of freedom: 2930
-DEAL:cg::Starting value 47.1699
-DEAL:cg::Convergence step 7 value 1.64836e-06
-DEAL::Testing FE_Q<3>(2)
-DEAL::Number of degrees of freedom: 1106
-DEAL:cg::Starting value 30.8707
-DEAL:cg::Convergence step 5 value 5.33421e-08
-DEAL::Testing FE_Q<3>(2)
-DEAL::Number of degrees of freedom: 4268
-DEAL:cg::Starting value 57.4891
-DEAL:cg::Convergence step 7 value 1.03521e-06
+DEAL:nothread::Testing FE_Q<2>(1)
+DEAL:nothread::Number of degrees of freedom: 507
+DEAL:nothread:cg::Starting value 21.9317
+DEAL:nothread:cg::Convergence step 5 value 7.81205e-07
+DEAL:threaded::Testing FE_Q<2>(1)
+DEAL:threaded::Number of degrees of freedom: 507
+DEAL:threaded:cg::Starting value 21.9317
+DEAL:threaded:cg::Convergence step 5 value 7.81205e-07
+DEAL:nothread::Testing FE_Q<2>(1)
+DEAL:nothread::Number of degrees of freedom: 881
+DEAL:nothread:cg::Starting value 27.7669
+DEAL:nothread:cg::Convergence step 6 value 9.24482e-07
+DEAL:threaded::Testing FE_Q<2>(1)
+DEAL:threaded::Number of degrees of freedom: 881
+DEAL:threaded:cg::Starting value 27.7669
+DEAL:threaded:cg::Convergence step 6 value 9.24482e-07
+DEAL:nothread::Testing FE_Q<2>(1)
+DEAL:nothread::Number of degrees of freedom: 2147
+DEAL:nothread:cg::Starting value 43.2551
+DEAL:nothread:cg::Convergence step 6 value 1.14520e-06
+DEAL:threaded::Testing FE_Q<2>(1)
+DEAL:threaded::Number of degrees of freedom: 2147
+DEAL:threaded:cg::Starting value 43.2551
+DEAL:threaded:cg::Convergence step 6 value 1.14520e-06
+DEAL:nothread::Testing FE_Q<2>(1)
+DEAL:nothread::Number of degrees of freedom: 6517
+DEAL:nothread:cg::Starting value 76.9220
+DEAL:nothread:cg::Convergence step 6 value 1.55925e-06
+DEAL:threaded::Testing FE_Q<2>(1)
+DEAL:threaded::Number of degrees of freedom: 6517
+DEAL:threaded:cg::Starting value 76.9220
+DEAL:threaded:cg::Convergence step 6 value 1.55925e-06
+DEAL:nothread::Testing FE_Q<2>(3)
+DEAL:nothread::Number of degrees of freedom: 4259
+DEAL:nothread:cg::Starting value 64.2573
+DEAL:nothread:cg::Convergence step 5 value 1.63125e-06
+DEAL:threaded::Testing FE_Q<2>(3)
+DEAL:threaded::Number of degrees of freedom: 4259
+DEAL:threaded:cg::Starting value 64.2573
+DEAL:threaded:cg::Convergence step 5 value 1.63120e-06
+DEAL:nothread::Testing FE_Q<2>(3)
+DEAL:nothread::Number of degrees of freedom: 7457
+DEAL:nothread:cg::Starting value 83.1084
+DEAL:nothread:cg::Convergence step 5 value 7.69258e-06
+DEAL:threaded::Testing FE_Q<2>(3)
+DEAL:threaded::Number of degrees of freedom: 7457
+DEAL:threaded:cg::Starting value 83.1084
+DEAL:threaded:cg::Convergence step 5 value 7.69315e-06
+DEAL:nothread::Testing FE_Q<2>(3)
+DEAL:nothread::Number of degrees of freedom: 18535
+DEAL:nothread:cg::Starting value 130.977
+DEAL:nothread:cg::Convergence step 5 value 6.04397e-06
+DEAL:threaded::Testing FE_Q<2>(3)
+DEAL:threaded::Number of degrees of freedom: 18535
+DEAL:threaded:cg::Starting value 130.977
+DEAL:threaded:cg::Convergence step 5 value 6.04357e-06
+DEAL:nothread::Testing FE_Q<3>(1)
+DEAL:nothread::Number of degrees of freedom: 186
+DEAL:nothread:cg::Starting value 12.3693
+DEAL:nothread:cg::Convergence step 3 value 8.40297e-08
+DEAL:threaded::Testing FE_Q<3>(1)
+DEAL:threaded::Number of degrees of freedom: 186
+DEAL:threaded:cg::Starting value 12.3693
+DEAL:threaded:cg::Convergence step 3 value 8.40297e-08
+DEAL:nothread::Testing FE_Q<3>(1)
+DEAL:nothread::Number of degrees of freedom: 648
+DEAL:nothread:cg::Starting value 21.1424
+DEAL:nothread:cg::Convergence step 6 value 1.05183e-07
+DEAL:threaded::Testing FE_Q<3>(1)
+DEAL:threaded::Number of degrees of freedom: 648
+DEAL:threaded:cg::Starting value 21.1424
+DEAL:threaded:cg::Convergence step 6 value 1.05183e-07
+DEAL:nothread::Testing FE_Q<3>(1)
+DEAL:nothread::Number of degrees of freedom: 2930
+DEAL:nothread:cg::Starting value 47.1699
+DEAL:nothread:cg::Convergence step 7 value 1.64836e-06
+DEAL:threaded::Testing FE_Q<3>(1)
+DEAL:threaded::Number of degrees of freedom: 2930
+DEAL:threaded:cg::Starting value 47.1699
+DEAL:threaded:cg::Convergence step 7 value 1.64836e-06
+DEAL:nothread::Testing FE_Q<3>(2)
+DEAL:nothread::Number of degrees of freedom: 1106
+DEAL:nothread:cg::Starting value 30.8707
+DEAL:nothread:cg::Convergence step 5 value 5.33407e-08
+DEAL:threaded::Testing FE_Q<3>(2)
+DEAL:threaded::Number of degrees of freedom: 1106
+DEAL:threaded:cg::Starting value 30.8707
+DEAL:threaded:cg::Convergence step 5 value 5.33448e-08
+DEAL:nothread::Testing FE_Q<3>(2)
+DEAL:nothread::Number of degrees of freedom: 4268
+DEAL:nothread:cg::Starting value 57.4891
+DEAL:nothread:cg::Convergence step 7 value 1.03521e-06
+DEAL:threaded::Testing FE_Q<3>(2)
+DEAL:threaded::Number of degrees of freedom: 4268
+DEAL:threaded:cg::Starting value 57.4891
+DEAL:threaded:cg::Convergence step 7 value 1.03521e-06