From d1d928c685837bf52e5718fe3483cd5674c2b59b Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Tue, 20 Nov 2012 15:31:31 +0000 Subject: [PATCH] Test the PETScWrappers::SparseDirectMumps interface. git-svn-id: https://svn.dealii.org/trunk@27593 0785d39b-7218-0410-832d-ea1e28bc413d --- tests/petsc/sparse_direct_mumps.cc | 94 +++++++++++++++++++++ tests/petsc/sparse_direct_mumps/cmp/generic | 4 + 2 files changed, 98 insertions(+) create mode 100644 tests/petsc/sparse_direct_mumps.cc create mode 100644 tests/petsc/sparse_direct_mumps/cmp/generic diff --git a/tests/petsc/sparse_direct_mumps.cc b/tests/petsc/sparse_direct_mumps.cc new file mode 100644 index 0000000000..c368f0eb95 --- /dev/null +++ b/tests/petsc/sparse_direct_mumps.cc @@ -0,0 +1,94 @@ +//---------------------------- sparse_direct_mumps.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. +// +//---------------------------- sparse_direct_mumps.cc --------------------------- + +// test the PETSc SparseDirectMumps solver + + +#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; +} + + +int main(int argc, char **argv) +{ + std::ofstream logfile("sparse_direct_mumps/output"); + deallog.attach(logfile); + deallog << std::setprecision(4); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + PetscInitialize(&argc,&argv,0,0); + { + 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); + u = 0.; + f = 1.; + A.compress (); + f.compress (); + u.compress (); + + SolverControl cn; + PETScWrappers::SparseDirectMUMPS solver(cn); +// solver.set_symmetric_mode(true); + solver.solve(A,u,f); + + PETScWrappers::Vector tmp(dim); + deallog << "residual = " << A.residual (tmp, u, f) + << std::endl; + } + PetscFinalize (); +} + diff --git a/tests/petsc/sparse_direct_mumps/cmp/generic b/tests/petsc/sparse_direct_mumps/cmp/generic new file mode 100644 index 0000000000..932860edcc --- /dev/null +++ b/tests/petsc/sparse_direct_mumps/cmp/generic @@ -0,0 +1,4 @@ + +DEAL::Size 32 Unknowns 961 +DEAL::Convergence step 1 value 0 +DEAL::residual = 0 -- 2.39.5