From: Wolfgang Bangerth Date: Tue, 20 Nov 2012 02:47:04 +0000 (+0000) Subject: Add another test. X-Git-Tag: v8.0.0~1811 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=c875f314ab0cdbf681a06377060535cc387a1c7d;p=dealii.git Add another test. git-svn-id: https://svn.dealii.org/trunk@27579 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/tests/petsc/solver_03_precondition_lu.cc b/tests/petsc/solver_03_precondition_lu.cc new file mode 100644 index 0000000000..bbcd02309c --- /dev/null +++ b/tests/petsc/solver_03_precondition_lu.cc @@ -0,0 +1,97 @@ +//---------------------------- petsc_solver_03_precondition_lu.cc --------------------------- +// $Id$ +// Version: $Name$ +// +// Copyright (C) 2004, 2005, 2010, 2012 by the deal.II authors +// +// This file is subject to QPL and may not be distributed +// without copyright and license information. Please refer +// to the file deal.II/doc/license.html for the text and +// further information on this license. +// +//---------------------------- petsc_solver_03_precondition_lu.cc --------------------------- + +// test the PETSc CG solver with the PreconditionLU +// preconditioner. This should converge in exactly one iteration + + +#include "../tests.h" +#include "../lac/testmatrix.h" +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +template +void +check_solve( SOLVER& solver, const MATRIX& A, + VECTOR& u, VECTOR& f, const PRECONDITION& P) +{ + deallog << "Solver type: " << typeid(solver).name() << std::endl; + + u = 0.; + f = 1.; + try + { + solver.solve(A,u,f,P); + } + catch (std::exception& e) + { + deallog << e.what() << std::endl; + abort (); + } + + deallog << "Solver stopped after " << solver.control().last_step() + << " iterations" << std::endl; + + // we use a direct solver as a + // preconditioner. this should + // converge in one step! + Assert (solver.control().last_step() == 1, + ExcInternalError()); +} + + +int main(int argc, char **argv) +{ + std::ofstream logfile("solver_03_precondition_lu/output"); + deallog.attach(logfile); + deallog << std::setprecision(4); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + PetscInitialize(&argc,&argv,0,0); + { + SolverControl control(100, 1.e-3); + + const unsigned int size = 32; + unsigned int dim = (size-1)*(size-1); + + deallog << "Size " << size << " Unknowns " << dim << std::endl; + + // Make matrix + FDMatrix testproblem(size, size); + PETScWrappers::SparseMatrix A(dim, dim, 5); + testproblem.five_point(A); + + PETScWrappers::Vector f(dim); + PETScWrappers::Vector u(dim); + f = 1.; + A.compress (); + f.compress (); + u.compress (); + + PETScWrappers::SolverCG solver(control); + PETScWrappers::PreconditionLU preconditioner(A); + check_solve (solver, A,u,f, preconditioner); + } + PetscFinalize (); +} + diff --git a/tests/petsc/solver_03_precondition_lu/cmp/generic b/tests/petsc/solver_03_precondition_lu/cmp/generic new file mode 100644 index 0000000000..071b4be96e --- /dev/null +++ b/tests/petsc/solver_03_precondition_lu/cmp/generic @@ -0,0 +1,6 @@ + +DEAL::Size 32 Unknowns 961 +DEAL::Solver type: N6dealii13PETScWrappers8SolverCGE +DEAL::Starting value 1351. +DEAL::Convergence step 1 value 0 +DEAL::Solver stopped after 1 iterations