From: David Wells Date: Mon, 7 Jun 2021 18:06:54 +0000 (-0400) Subject: Use initlog in some tests. X-Git-Tag: v9.4.0-rc1~1254^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=5ccbb0cc74935d230f9b247a908e216f586e65a3;p=dealii.git Use initlog in some tests. We forgot to call mpi_initlog() in these tests, which caused the tests to compare stderr + stdout as output - this is why we had a bunch of problems with testers collecting additional MPI error output. --- diff --git a/tests/simplex/step-17.cc b/tests/simplex/step-17.cc index 0e26669268..fe0e208ce2 100644 --- a/tests/simplex/step-17.cc +++ b/tests/simplex/step-17.cc @@ -74,6 +74,8 @@ #include +#include "../tests.h" + // Flag to toggle between hexes and simplices. // #define HEX @@ -107,8 +109,6 @@ namespace Step17 const unsigned int n_mpi_processes; const unsigned int this_mpi_process; - ConditionalOStream pcout; - #ifdef HEX MappingQGeneric mapping; #else @@ -174,7 +174,6 @@ namespace Step17 : mpi_communicator(MPI_COMM_WORLD) , n_mpi_processes(Utilities::MPI::n_mpi_processes(mpi_communicator)) , this_mpi_process(Utilities::MPI::this_mpi_process(mpi_communicator)) - , pcout(std::cout, (this_mpi_process == 0)) #ifdef HEX , mapping(1) , fe(FE_Q(1), dim) @@ -329,9 +328,17 @@ namespace Step17 SolverControl solver_control(solution.size(), 1e-8 * system_rhs.l2_norm()); PETScWrappers::SolverCG cg(solver_control, mpi_communicator); + // PreconditionBlockJacobi depends on the processor count PETScWrappers::PreconditionBlockJacobi preconditioner(system_matrix); + const unsigned int lower = n_mpi_processes == 1 ? 8 : 26; + const unsigned int upper = n_mpi_processes == 1 ? 12 : 30; + + check_solver_within_range( + cg.solve(system_matrix, solution, system_rhs, preconditioner), + solver_control.last_step(), + lower, + upper); - cg.solve(system_matrix, solution, system_rhs, preconditioner); Vector localized_solution(solution); hanging_node_constraints.distribute(localized_solution); @@ -420,25 +427,21 @@ namespace Step17 triangulation, std::vector(dim, n_subdivisions), a, b); #endif - pcout << " Number of active cells: " - << triangulation.n_active_cells() << std::endl; + deallog << " Number of active cells: " + << triangulation.n_active_cells() << std::endl; setup_system(); - pcout << " Number of degrees of freedom: " << dof_handler.n_dofs() - << " (by partition:"; + deallog << " Number of degrees of freedom: " << dof_handler.n_dofs() + << " (by partition:"; for (unsigned int p = 0; p < n_mpi_processes; ++p) - pcout << (p == 0 ? ' ' : '+') - << (DoFTools::count_dofs_with_subdomain_association(dof_handler, - p)); - pcout << ")" << std::endl; + deallog << (p == 0 ? ' ' : '+') + << (DoFTools::count_dofs_with_subdomain_association(dof_handler, + p)); + deallog << ")" << std::endl; assemble_system(); - const unsigned int n_iterations = solve(); - - pcout << " Solver converged in " << n_iterations << " iterations." - << std::endl; - + solve(); output_results(); } } // namespace Step17 @@ -447,13 +450,14 @@ namespace Step17 int main(int argc, char **argv) { + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + mpi_initlog(); + try { using namespace dealii; using namespace Step17; - Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); - ElasticProblem<2> elastic_problem; elastic_problem.run(); } diff --git a/tests/simplex/step-17.mpirun=1.with_petsc=true.output b/tests/simplex/step-17.mpirun=1.with_petsc=true.output index 0a47900732..7ef60db63d 100644 --- a/tests/simplex/step-17.mpirun=1.with_petsc=true.output +++ b/tests/simplex/step-17.mpirun=1.with_petsc=true.output @@ -1,3 +1,4 @@ - Number of active cells: 200 - Number of degrees of freedom: 242 (by partition: 242) - Solver converged in 10 iterations. + +DEAL:: Number of active cells: 200 +DEAL:: Number of degrees of freedom: 242 (by partition: 242) +DEAL::Solver stopped within 8 - 12 iterations diff --git a/tests/simplex/step-17.mpirun=4.with_petsc=true.output b/tests/simplex/step-17.mpirun=4.with_petsc=true.output index 956179849e..49b50d8304 100644 --- a/tests/simplex/step-17.mpirun=4.with_petsc=true.output +++ b/tests/simplex/step-17.mpirun=4.with_petsc=true.output @@ -1,3 +1,4 @@ - Number of active cells: 200 - Number of degrees of freedom: 242 (by partition: 74+58+58+52) - Solver converged in 28 iterations. + +DEAL:: Number of active cells: 200 +DEAL:: Number of degrees of freedom: 242 (by partition: 74+58+58+52) +DEAL::Solver stopped within 26 - 30 iterations diff --git a/tests/simplex/step-40.cc b/tests/simplex/step-40.cc index 54c5fc2c2c..7a670fbe41 100644 --- a/tests/simplex/step-40.cc +++ b/tests/simplex/step-40.cc @@ -85,6 +85,8 @@ namespace LA #include #include +#include "../tests.h" + // Flag to toggle between hex and simplex. // #define HEX @@ -136,8 +138,6 @@ namespace Step40 LA::MPI::SparseMatrix system_matrix; LA::MPI::Vector locally_relevant_solution; LA::MPI::Vector system_rhs; - - ConditionalOStream pcout; }; template @@ -155,8 +155,6 @@ namespace Step40 #endif , fe(2) , dof_handler(triangulation) - , pcout(std::cout, - (Utilities::MPI::this_mpi_process(mpi_communicator) == 0)) {} template @@ -287,10 +285,13 @@ namespace Step40 #endif preconditioner.initialize(system_matrix, data); - solver.solve(system_matrix, - completely_distributed_solution, - system_rhs, - preconditioner); + check_solver_within_range(solver.solve(system_matrix, + completely_distributed_solution, + system_rhs, + preconditioner), + solver_control.last_step(), + 5, + 9); constraints.distribute(completely_distributed_solution); @@ -312,10 +313,10 @@ namespace Step40 quadrature_formula, VectorTools::NormType::L2_norm); - pcout << VectorTools::compute_global_error(triangulation, - difference, - VectorTools::NormType::L2_norm) - << std::endl; + deallog << VectorTools::compute_global_error(triangulation, + difference, + VectorTools::NormType::L2_norm) + << std::endl; } template @@ -341,14 +342,14 @@ namespace Step40 void LaplaceProblem::run() { - pcout << "Running with " + deallog << "Running with " #ifdef USE_PETSC_LA - << "PETSc" + << "PETSc" #else - << "Trilinos" + << "Trilinos" #endif - << " on " << Utilities::MPI::n_mpi_processes(mpi_communicator) - << " MPI rank(s)..." << std::endl; + << " on " << Utilities::MPI::n_mpi_processes(mpi_communicator) + << " MPI rank(s)..." << std::endl; const unsigned int n_subdivisions = 16; Point a, b; @@ -394,10 +395,10 @@ namespace Step40 setup_system(); - pcout << " Number of active cells: " - << triangulation.n_global_active_cells() << std::endl - << " Number of degrees of freedom: " << dof_handler.n_dofs() - << std::endl; + deallog << " Number of active cells: " + << triangulation.n_global_active_cells() << std::endl + << " Number of degrees of freedom: " << dof_handler.n_dofs() + << std::endl; assemble_system(); solve(); @@ -413,12 +414,14 @@ namespace Step40 int main(int argc, char *argv[]) { + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + mpi_initlog(); + try { using namespace dealii; using namespace Step40; - Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); LaplaceProblem<2> laplace_problem_2d; laplace_problem_2d.run(); diff --git a/tests/simplex/step-40.mpirun=1.with_petsc=true.output b/tests/simplex/step-40.mpirun=1.with_petsc=true.output index ca81c16970..5693999f84 100644 --- a/tests/simplex/step-40.mpirun=1.with_petsc=true.output +++ b/tests/simplex/step-40.mpirun=1.with_petsc=true.output @@ -1,4 +1,6 @@ -Running with PETSc on 1 MPI rank(s)... - Number of active cells: 512 - Number of degrees of freedom: 1089 -0.0124965 + +DEAL::Running with PETSc on 1 MPI rank(s)... +DEAL:: Number of active cells: 512 +DEAL:: Number of degrees of freedom: 1089 +DEAL::Solver stopped within 5 - 9 iterations +DEAL::0.0124965 diff --git a/tests/simplex/step-40.mpirun=4.with_petsc=true.with_metis=true.output b/tests/simplex/step-40.mpirun=4.with_petsc=true.with_metis=true.output index d476df598c..bff795f87f 100644 --- a/tests/simplex/step-40.mpirun=4.with_petsc=true.with_metis=true.output +++ b/tests/simplex/step-40.mpirun=4.with_petsc=true.with_metis=true.output @@ -1,4 +1,6 @@ -Running with PETSc on 4 MPI rank(s)... - Number of active cells: 512 - Number of degrees of freedom: 1089 -0.0124965 + +DEAL::Running with PETSc on 4 MPI rank(s)... +DEAL:: Number of active cells: 512 +DEAL:: Number of degrees of freedom: 1089 +DEAL::Solver stopped within 5 - 9 iterations +DEAL::0.0124965