From: David Wells Date: Mon, 17 Feb 2025 16:22:50 +0000 (-0500) Subject: Fix the frequently-failing trilinos/solver_control_06 test. X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=be28dcd04008880cc5a2977b9a6b92af6d066f4b;p=dealii.git Fix the frequently-failing trilinos/solver_control_06 test. The main problem here is // Trilinos dumps the output into std::cout // We catch this output and it is written to the stdout logfile // Since we're interested in this output we read it back in and // write parts of it to the logstream // We can only do this reliably if we are not running in parallel (sometimes // stdout is not written yet otherwise) if (Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) == 1) { std::ifstream inputfile; inputfile.open("stdout"); Assert(inputfile.good() && inputfile.is_open(), ExcIO()); std::string line; const std::string key = "*****"; while (std::getline(inputfile, line)) { if (line.find(key) != std::string::npos) deallog << line << std::endl; } inputfile.close(); } deallog << "OK" << std::endl; which falsely assumes that, at the point we hit that check, all buffered output has been written. In general we cannot guarantee this until the end of the program since we are actively, during program execution, piping stdout and sterr to output files. We don't really need to check this output (Trilinos' residuals) anyway - we already know, if the solver converges, that we got close enough. Hence, this patch replaces all of this with our standard check_solver_within_range() macro which is vastly more robust. --- diff --git a/tests/trilinos/solver_control_06.cc b/tests/trilinos/solver_control_06.cc index d06fc9b425..572869dcac 100644 --- a/tests/trilinos/solver_control_06.cc +++ b/tests/trilinos/solver_control_06.cc @@ -77,17 +77,17 @@ private: assemble_system(); void - solve_base(); + solve_base(const unsigned int cycle); void - solve_cg(); + solve_cg(const unsigned int cycle); void - solve_cgs(); + solve_cgs(const unsigned int cycle); void - solve_gmres(); + solve_gmres(const unsigned int cycle); void - solve_bicgstab(); + solve_bicgstab(const unsigned int cycle); void - solve_tfqmr(); + solve_tfqmr(const unsigned int cycle); void refine_grid(); @@ -113,8 +113,9 @@ private: LA::MPI::Vector locally_relevant_solution; LA::MPI::Vector system_rhs; - ConditionalOStream pcout; - TimerOutput timer; + ConditionalOStream pcout_timer; + + TimerOutput timer; }; Test_Solver_Output::Test_Solver_Output() @@ -127,11 +128,8 @@ Test_Solver_Output::Test_Solver_Output() Triangulation<2>::smoothing_on_coarsening)) , dof_handler(triangulation) , fe(1) - , pcout(std::cout, (Utilities::MPI::this_mpi_process(mpi_comm) == 0)) - , - // pcout(deallog.get_file_stream(), - // (Utilities::MPI::this_mpi_process(mpi_comm) == 0)), - timer(mpi_comm, pcout, TimerOutput::never, TimerOutput::wall_times) + , pcout_timer(std::cout, (Utilities::MPI::this_mpi_process(mpi_comm) == 0)) + , timer(mpi_comm, pcout_timer, TimerOutput::never, TimerOutput::wall_times) {} Test_Solver_Output::~Test_Solver_Output() @@ -145,7 +143,7 @@ Test_Solver_Output::run() const unsigned int n_cycles = 2; for (unsigned int cycle = 0; cycle < n_cycles; ++cycle) { - pcout << " Cycle: " << cycle << std::endl; + deallog << " Cycle: " << cycle << std::endl; if (cycle == 0) { make_grid(); @@ -157,18 +155,18 @@ Test_Solver_Output::run() 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_base(); - solve_cg(); - solve_cgs(); - solve_gmres(); - solve_bicgstab(); - solve_tfqmr(); + solve_base(cycle); + solve_cg(cycle); + solve_cgs(cycle); + solve_gmres(cycle); + solve_bicgstab(cycle); + solve_tfqmr(cycle); // { // TimerOutput::Scope t(timer, "output"); @@ -178,7 +176,7 @@ Test_Solver_Output::run() // timer.print_summary(); // timer.reset(); - pcout << std::endl << std::endl << std::endl; + deallog << std::endl << std::endl << std::endl; } } @@ -286,10 +284,10 @@ Test_Solver_Output::assemble_system() } void -Test_Solver_Output::solve_base() +Test_Solver_Output::solve_base(const unsigned int cycle) { TimerOutput::Scope t(timer, "solve_base"); - pcout << "Solving using SolverBase" << std::endl; + deallog << "Solving using SolverBase" << std::endl; LA::MPI::Vector completely_distributed_solution(locally_owned_dofs, mpi_comm); @@ -302,23 +300,23 @@ Test_Solver_Output::solve_base() TrilinosWrappers::SolverBase solver(TrilinosWrappers::SolverBase::cg, solver_control, solver_data); - solver.solve(system_matrix, - completely_distributed_solution, - system_rhs, - prec); - - pcout << " Solved in " << solver_control.last_step() << " iterations." - << std::endl; + check_solver_within_range(solver.solve(system_matrix, + completely_distributed_solution, + system_rhs, + prec), + solver_control.last_step(), + (cycle == 0 ? 1 : 11), + (cycle == 0 ? 2 : 15)); constraints.distribute(completely_distributed_solution); locally_relevant_solution = completely_distributed_solution; } void -Test_Solver_Output::solve_cg() +Test_Solver_Output::solve_cg(const unsigned int cycle) { TimerOutput::Scope t(timer, "solve_cg"); - pcout << "Solving using SolverCG" << std::endl; + deallog << "Solving using SolverCG" << std::endl; LA::MPI::Vector completely_distributed_solution(locally_owned_dofs, mpi_comm); @@ -329,23 +327,23 @@ Test_Solver_Output::solve_cg() SolverControl solver_control(100, 1e-12); LA::SolverCG::AdditionalData solver_data(true); LA::SolverCG solver(solver_control, solver_data); - solver.solve(system_matrix, - completely_distributed_solution, - system_rhs, - prec); - - pcout << " Solved in " << solver_control.last_step() << " iterations." - << std::endl; + check_solver_within_range(solver.solve(system_matrix, + completely_distributed_solution, + system_rhs, + prec), + solver_control.last_step(), + (cycle == 0 ? 1 : 11), + (cycle == 0 ? 2 : 15)); constraints.distribute(completely_distributed_solution); locally_relevant_solution = completely_distributed_solution; } void -Test_Solver_Output::solve_cgs() +Test_Solver_Output::solve_cgs(const unsigned int cycle) { TimerOutput::Scope t(timer, "solve_cgs"); - pcout << "Solving using SolverCGS" << std::endl; + deallog << "Solving using SolverCGS" << std::endl; LA::MPI::Vector completely_distributed_solution(locally_owned_dofs, mpi_comm); @@ -356,23 +354,23 @@ Test_Solver_Output::solve_cgs() SolverControl solver_control(100, 1e-12); TrilinosWrappers::SolverCGS::AdditionalData solver_data(true); TrilinosWrappers::SolverCGS solver(solver_control, solver_data); - solver.solve(system_matrix, - completely_distributed_solution, - system_rhs, - prec); - - pcout << " Solved in " << solver_control.last_step() << " iterations." - << std::endl; + check_solver_within_range(solver.solve(system_matrix, + completely_distributed_solution, + system_rhs, + prec), + solver_control.last_step(), + (cycle == 0 ? 1 : 5), + (cycle == 0 ? 2 : 9)); constraints.distribute(completely_distributed_solution); locally_relevant_solution = completely_distributed_solution; } void -Test_Solver_Output::solve_gmres() +Test_Solver_Output::solve_gmres(const unsigned int cycle) { TimerOutput::Scope t(timer, "solve_gmres"); - pcout << "Solving using SolverGMRES" << std::endl; + deallog << "Solving using SolverGMRES" << std::endl; LA::MPI::Vector completely_distributed_solution(locally_owned_dofs, mpi_comm); @@ -383,23 +381,23 @@ Test_Solver_Output::solve_gmres() SolverControl solver_control(100, 1e-12); LA::SolverGMRES::AdditionalData solver_data(true, 25); LA::SolverGMRES solver(solver_control, solver_data); - solver.solve(system_matrix, - completely_distributed_solution, - system_rhs, - prec); - - pcout << " Solved in " << solver_control.last_step() << " iterations." - << std::endl; + check_solver_within_range(solver.solve(system_matrix, + completely_distributed_solution, + system_rhs, + prec), + solver_control.last_step(), + (cycle == 0 ? 1 : 11), + (cycle == 0 ? 2 : 15)); constraints.distribute(completely_distributed_solution); locally_relevant_solution = completely_distributed_solution; } void -Test_Solver_Output::solve_bicgstab() +Test_Solver_Output::solve_bicgstab(const unsigned int cycle) { TimerOutput::Scope t(timer, "solve_bicgstab"); - pcout << "Solving using SolverBicgstab" << std::endl; + deallog << "Solving using SolverBicgstab" << std::endl; LA::MPI::Vector completely_distributed_solution(locally_owned_dofs, mpi_comm); @@ -410,23 +408,23 @@ Test_Solver_Output::solve_bicgstab() SolverControl solver_control(100, 1e-12); TrilinosWrappers::SolverBicgstab::AdditionalData solver_data(true); TrilinosWrappers::SolverBicgstab solver(solver_control, solver_data); - solver.solve(system_matrix, - completely_distributed_solution, - system_rhs, - prec); - - pcout << " Solved in " << solver_control.last_step() << " iterations." - << std::endl; + check_solver_within_range(solver.solve(system_matrix, + completely_distributed_solution, + system_rhs, + prec), + solver_control.last_step(), + (cycle == 0 ? 1 : 5), + (cycle == 0 ? 2 : 9)); constraints.distribute(completely_distributed_solution); locally_relevant_solution = completely_distributed_solution; } void -Test_Solver_Output::solve_tfqmr() +Test_Solver_Output::solve_tfqmr(const unsigned int cycle) { TimerOutput::Scope t(timer, "solve_tfqmr"); - pcout << "Solving using SolverTFQMR" << std::endl; + deallog << "Solving using SolverTFQMR" << std::endl; LA::MPI::Vector completely_distributed_solution(locally_owned_dofs, mpi_comm); @@ -437,13 +435,13 @@ Test_Solver_Output::solve_tfqmr() SolverControl solver_control(100, 1e-12); TrilinosWrappers::SolverTFQMR::AdditionalData solver_data(true); TrilinosWrappers::SolverTFQMR solver(solver_control, solver_data); - solver.solve(system_matrix, - completely_distributed_solution, - system_rhs, - prec); - - pcout << " Solved in " << solver_control.last_step() << " iterations." - << std::endl; + check_solver_within_range(solver.solve(system_matrix, + completely_distributed_solution, + system_rhs, + prec), + solver_control.last_step(), + (cycle == 0 ? 1 : 5), + (cycle == 0 ? 2 : 9)); constraints.distribute(completely_distributed_solution); locally_relevant_solution = completely_distributed_solution; @@ -510,28 +508,4 @@ main(int argc, char *argv[]) Test_Solver_Output problem; problem.run(); - - // Trilinos dumps the output into std::cout - // We catch this output and it is written to the stdout logfile - // Since we're interested in this output we read it back in and - // write parts of it to the logstream - // We can only do this reliably if we are not running in parallel (sometimes - // stdout is not written yet otherwise) - if (Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) == 1) - { - std::ifstream inputfile; - inputfile.open("stdout"); - Assert(inputfile.good() && inputfile.is_open(), ExcIO()); - std::string line; - const std::string key = "*****"; - while (std::getline(inputfile, line)) - { - if (line.find(key) != std::string::npos) - deallog << line << std::endl; - } - inputfile.close(); - } - deallog << "OK" << std::endl; - - return 0; } diff --git a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output index 9a060367ff..255595f6b3 100644 --- a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output +++ b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output @@ -1,74 +1,37 @@ -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL::OK +DEAL:: Cycle: 0 +DEAL:: Number of active cells: 1024 +DEAL:: Number of degrees of freedom: 1089 +DEAL::Solving using SolverBase +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverCG +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverCGS +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverGMRES +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverBicgstab +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverTFQMR +DEAL::Solver stopped within 1 - 2 iterations +DEAL:: +DEAL:: +DEAL:: +DEAL:: Cycle: 1 +DEAL:: Number of active cells: 1960 +DEAL:: Number of degrees of freedom: 2153 +DEAL::Solving using SolverBase +DEAL::Solver stopped within 11 - 15 iterations +DEAL::Solving using SolverCG +DEAL::Solver stopped within 11 - 15 iterations +DEAL::Solving using SolverCGS +DEAL::Solver stopped within 5 - 9 iterations +DEAL::Solving using SolverGMRES +DEAL::Solver stopped within 11 - 15 iterations +DEAL::Solving using SolverBicgstab +DEAL::Solver stopped within 5 - 9 iterations +DEAL::Solving using SolverTFQMR +DEAL::Solver stopped within 5 - 9 iterations +DEAL:: +DEAL:: +DEAL:: diff --git a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output.2 b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output.2 deleted file mode 100644 index 0fd8fc12f0..0000000000 --- a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output.2 +++ /dev/null @@ -1,2 +0,0 @@ - -DEAL::OK diff --git a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output.3 b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output.3 deleted file mode 100644 index 85b83dd2ff..0000000000 --- a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=1.output.3 +++ /dev/null @@ -1,74 +0,0 @@ - -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL::OK diff --git a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output index 9a060367ff..255595f6b3 100644 --- a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output +++ b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output @@ -1,74 +1,37 @@ -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=1, ~/Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=2, Cheby_pre0/Cheby_post0, ~/Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL::OK +DEAL:: Cycle: 0 +DEAL:: Number of active cells: 1024 +DEAL:: Number of degrees of freedom: 1089 +DEAL::Solving using SolverBase +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverCG +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverCGS +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverGMRES +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverBicgstab +DEAL::Solver stopped within 1 - 2 iterations +DEAL::Solving using SolverTFQMR +DEAL::Solver stopped within 1 - 2 iterations +DEAL:: +DEAL:: +DEAL:: +DEAL:: Cycle: 1 +DEAL:: Number of active cells: 1960 +DEAL:: Number of degrees of freedom: 2153 +DEAL::Solving using SolverBase +DEAL::Solver stopped within 11 - 15 iterations +DEAL::Solving using SolverCG +DEAL::Solver stopped within 11 - 15 iterations +DEAL::Solving using SolverCGS +DEAL::Solver stopped within 5 - 9 iterations +DEAL::Solving using SolverGMRES +DEAL::Solver stopped within 11 - 15 iterations +DEAL::Solving using SolverBicgstab +DEAL::Solver stopped within 5 - 9 iterations +DEAL::Solving using SolverTFQMR +DEAL::Solver stopped within 5 - 9 iterations +DEAL:: +DEAL:: +DEAL:: diff --git a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output.2 b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output.2 deleted file mode 100644 index 0fd8fc12f0..0000000000 --- a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output.2 +++ /dev/null @@ -1,2 +0,0 @@ - -DEAL::OK diff --git a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output.3 b/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output.3 deleted file mode 100644 index 85b83dd2ff..0000000000 --- a/tests/trilinos/solver_control_06.with_p4est=true.mpirun=2.output.3 +++ /dev/null @@ -1,74 +0,0 @@ - -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=1, /Amesos_KLU_0, ) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CG solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned CGS solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned GMRES solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned BICGSTAB solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL:: ******************************************************* -DEAL:: ***** Problem: Epetra::CrsMatrix -DEAL:: ***** Preconditioned TFQMR solution -DEAL:: ***** ML (L=2, /Cheby_post0, /Amesos_KLU_1) -DEAL:: ***** No scaling -DEAL:: ******************************************************* -DEAL::OK